Panzer  Version of the Day
Panzer_L2Projection_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // @HEADER
3 
4 #ifndef PANZER_L2_PROJECTION_IMPL_HPP
5 #define PANZER_L2_PROJECTION_IMPL_HPP
6 
8 #include "Tpetra_CrsGraph.hpp"
9 #include "Tpetra_MultiVector.hpp"
10 #include "Tpetra_CrsMatrix.hpp"
13 #include "Panzer_TpetraLinearObjFactory.hpp"
16 #include "Panzer_DOFManagerFactory.hpp"
17 #include "Panzer_BlockedDOFManagerFactory.hpp"
21 #include "Panzer_Workset.hpp"
22 
23 namespace panzer {
24 
25  template<typename LO, typename GO>
27  const panzer::IntegrationDescriptor& integrationDescriptor,
28  const Teuchos::RCP<const Teuchos::MpiComm<int>>& comm,
29  const Teuchos::RCP<const panzer::ConnManager<LO,GO>>& connManager,
30  const std::vector<std::string>& elementBlockNames,
31  const Teuchos::RCP<panzer::WorksetContainer> worksetContainer)
32  {
33  targetBasisDescriptor_ = targetBasis;
34  integrationDescriptor_ = integrationDescriptor;
35  comm_ = comm;
36  connManager_ = connManager;
37  elementBlockNames_ = elementBlockNames;
38  worksetContainer_ = worksetContainer;
39  setupCalled_ = true;
40 
41  // Build target DOF Manager
42  targetGlobalIndexer_ =
43  Teuchos::rcp(new panzer::DOFManager<LO,GO>(Teuchos::rcp_const_cast<panzer::ConnManager<LO,GO>>(connManager),*(comm->getRawMpiComm())));
44 
45  // For hybrid mesh, blocks could have different topology
46  for (const auto& eBlock : elementBlockNames_) {
47  std::vector<shards::CellTopology> topologies;
48  connManager_->getElementBlockTopologies(topologies);
49  std::vector<std::string> ebNames;
50  connManager_->getElementBlockIds(ebNames);
51  const auto search = std::find(ebNames.cbegin(),ebNames.cend(),eBlock);
52  TEUCHOS_ASSERT(search != ebNames.cend());
53  const int index = std::distance(ebNames.cbegin(),search);
54  const auto& cellTopology = topologies[index];
55 
56  auto intrepidBasis = panzer::createIntrepid2Basis<PHX::Device,double,double>(targetBasisDescriptor_.getType(),
57  targetBasisDescriptor_.getOrder(),
58  cellTopology);
60  // Field name is the basis type
61  targetGlobalIndexer_->addField(eBlock,targetBasisDescriptor_.getType(),fieldPattern);
62  }
63 
64  targetGlobalIndexer_->buildGlobalUnknowns();
65 
66  // Check workset needs are correct
67  }
68 
69  template<typename LO, typename GO>
72  {return targetGlobalIndexer_;}
73 
74  template<typename LO, typename GO>
77  {
78  TEUCHOS_ASSERT(setupCalled_);
79 
80  // Allocate the owned matrix
81  std::vector<Teuchos::RCP<const panzer::UniqueGlobalIndexer<LO,GO>>> indexers;
82  indexers.push_back(targetGlobalIndexer_);
83 
85 
86  auto ownedMatrix = factory.getTpetraMatrix(0,0);
87  auto ghostedMatrix = factory.getGhostedTpetraMatrix(0,0);
88  ownedMatrix->resumeFill();
89  ownedMatrix->setAllToScalar(0.0);
90  ghostedMatrix->resumeFill();
91  ghostedMatrix->setAllToScalar(0.0);
92  PHX::Device::fence();
93 
94  auto M = ghostedMatrix->getLocalMatrix();
95  const int fieldIndex = targetGlobalIndexer_->getFieldNum(targetBasisDescriptor_.getType());
96 
97  const bool is_scalar = targetBasisDescriptor_.getType()=="HGrad" || targetBasisDescriptor_.getType()=="Const";
98 
99  // Loop over element blocks and fill mass matrix
100  if(is_scalar){
101  for (const auto& block : elementBlockNames_) {
102 
103  // Based on descriptor, currently assumes there should only be one workset
105  const auto& worksets = worksetContainer_->getWorksets(wd);
106 
107  for (const auto& workset : *worksets) {
108 
109  const auto& basisValues = workset.getBasisValues(targetBasisDescriptor_,integrationDescriptor_);
110 
111  const auto& unweightedBasis = basisValues.basis_scalar;
112  const auto& weightedBasis = basisValues.weighted_basis_scalar;
113 
114  // Offsets (this assumes UVM, need to fix)
115  const std::vector<LO>& offsets = targetGlobalIndexer_->getGIDFieldOffsets(block,fieldIndex);
116  PHX::View<LO*> kOffsets("MassMatrix: Offsets",offsets.size());
117  for(const auto& i : offsets)
118  kOffsets(i) = offsets[i];
119 
120  PHX::Device::fence();
121 
122  // Local Ids
123  Kokkos::View<LO**,PHX::Device> localIds("MassMatrix: LocalIds", workset.numOwnedCells()+workset.numGhostCells()+workset.numVirtualCells(),
124  targetGlobalIndexer_->getElementBlockGIDCount(block));
125 
126  // Remove the ghosted cell ids or the call to getElementLocalIds will spill array bounds
127  const auto cellLocalIdsNoGhost = Kokkos::subview(workset.cell_local_ids_k,std::make_pair(0,workset.numOwnedCells()));
128 
129  targetGlobalIndexer_->getElementLIDs(cellLocalIdsNoGhost,localIds);
130 
131  const int numBasisPoints = static_cast<int>(weightedBasis.extent(1));
132  Kokkos::parallel_for(workset.numOwnedCells(),KOKKOS_LAMBDA (const int& cell) {
133 
134  LO cLIDs[256];
135  const int numIds = static_cast<int>(localIds.extent(1));
136  for(int i=0;i<numIds;++i)
137  cLIDs[i] = localIds(cell,i);
138 
139  double vals[256];
140  const int numQP = static_cast<int>(unweightedBasis.extent(2));
141 
142  for (int qp=0; qp < numQP; ++qp) {
143  for (int row=0; row < numBasisPoints; ++row) {
144  int offset = kOffsets(row);
145  LO lid = localIds(cell,offset);
146 
147  for (int col=0; col < numIds; ++col)
148  vals[col] = unweightedBasis(cell,row,qp) * weightedBasis(cell,col,qp);
149 
150  M.sumIntoValues(lid,cLIDs,numIds,vals,true,true);
151  }
152  }
153 
154  });
155  }
156  }
157  } else {
158  for (const auto& block : elementBlockNames_) {
159 
160  // Based on descriptor, currently assumes there should only be one workset
162  const auto& worksets = worksetContainer_->getWorksets(wd);
163 
164  for (const auto& workset : *worksets) {
165 
166  const auto& basisValues = workset.getBasisValues(targetBasisDescriptor_,integrationDescriptor_);
167 
168  const auto& unweightedBasis = basisValues.basis_vector;
169  const auto& weightedBasis = basisValues.weighted_basis_vector;
170 
171  // Offsets (this assumes UVM, need to fix)
172  const std::vector<LO>& offsets = targetGlobalIndexer_->getGIDFieldOffsets(block,fieldIndex);
173  PHX::View<LO*> kOffsets("MassMatrix: Offsets",offsets.size());
174  for(const auto& i : offsets)
175  kOffsets(i) = offsets[i];
176 
177  PHX::Device::fence();
178 
179  // Local Ids
180  Kokkos::View<LO**,PHX::Device> localIds("MassMatrix: LocalIds", workset.numOwnedCells()+workset.numGhostCells()+workset.numVirtualCells(),
181  targetGlobalIndexer_->getElementBlockGIDCount(block));
182 
183  // Remove the ghosted cell ids or the call to getElementLocalIds will spill array bounds
184  const auto cellLocalIdsNoGhost = Kokkos::subview(workset.cell_local_ids_k,std::make_pair(0,workset.numOwnedCells()));
185 
186  targetGlobalIndexer_->getElementLIDs(cellLocalIdsNoGhost,localIds);
187 
188  const int numBasisPoints = static_cast<int>(weightedBasis.extent(1));
189  Kokkos::parallel_for(workset.numOwnedCells(),KOKKOS_LAMBDA (const int& cell) {
190 
191  LO cLIDs[256];
192  const int numIds = static_cast<int>(localIds.extent(1));
193  for(int i=0;i<numIds;++i)
194  cLIDs[i] = localIds(cell,i);
195 
196  double vals[256];
197  const int numQP = static_cast<int>(unweightedBasis.extent(2));
198 
199  for (int qp=0; qp < numQP; ++qp) {
200  for (int row=0; row < numBasisPoints; ++row) {
201  int offset = kOffsets(row);
202  LO lid = localIds(cell,offset);
203 
204  for (int col=0; col < numIds; ++col){
205  vals[col] = 0.0;
206  for(int dim=0; dim < weightedBasis.extent(3); dim++)
207  vals[col] += unweightedBasis(cell,row,qp,dim) * weightedBasis(cell,col,qp,dim);
208  }
209 
210  M.sumIntoValues(lid,cLIDs,numIds,vals,true,true);
211  }
212  }
213 
214  });
215  }
216  }
217  }
218  PHX::exec_space::fence();
219 
220  ghostedMatrix->fillComplete();
221  const auto exporter = factory.getGhostedExport(0);
222  ownedMatrix->doExport(*ghostedMatrix, *exporter, Tpetra::ADD);
223  ownedMatrix->fillComplete();
224 
225  return ownedMatrix;
226  }
227 
228  template<typename LO, typename GO>
231  {
232  using Teuchos::rcp;
233  const auto massMatrix = this->buildMassMatrix();
234  const auto lumpedMassMatrix = rcp(new Tpetra::MultiVector<double,LO,GO,Kokkos::Compat::KokkosDeviceWrapperNode<PHX::Device>>(massMatrix->getDomainMap(),1,true));
235  const auto tmp = rcp(new Tpetra::MultiVector<double,LO,GO,Kokkos::Compat::KokkosDeviceWrapperNode<PHX::Device>>(massMatrix->getRangeMap(),1,false));
236  tmp->putScalar(1.0);
237  massMatrix->apply(*tmp,*lumpedMassMatrix);
238  lumpedMassMatrix->reciprocal(*lumpedMassMatrix);
239  return lumpedMassMatrix;
240  }
241 
242  template<typename LO, typename GO>
245  const Teuchos::RCP<const Tpetra::Map<LO,GO,Kokkos::Compat::KokkosDeviceWrapperNode<PHX::Device>>>& inputOwnedSourceMap,
246  const std::string& sourceFieldName,
247  const panzer::BasisDescriptor& sourceBasisDescriptor,
248  const int directionIndex)
249  {
250  // *******************
251  // Create Retangular matrix (both ghosted and owned).
252  // *******************
253  using Teuchos::RCP;
254  using Teuchos::rcp;
255  using MapType = Tpetra::Map<LO,GO,Kokkos::Compat::KokkosDeviceWrapperNode<PHX::Device>>;
256  using GraphType = Tpetra::CrsGraph<LO,GO,Kokkos::Compat::KokkosDeviceWrapperNode<PHX::Device>>;
257  using ExportType = Tpetra::Export<LO,GO,Kokkos::Compat::KokkosDeviceWrapperNode<PHX::Device>>;
258  using MatrixType = Tpetra::CrsMatrix<double,LO,GO,Kokkos::Compat::KokkosDeviceWrapperNode<PHX::Device>>;
259 
260  // *******************
261  // Build ghosted graph
262  // *******************
263 
264  RCP<MapType> ghostedTargetMap;
265  {
266  std::vector<GO> indices;
267  targetGlobalIndexer_->getOwnedAndGhostedIndices(indices);
268  ghostedTargetMap = rcp(new MapType(Teuchos::OrdinalTraits<GO>::invalid(),indices,0,comm_));
269  }
270 
271  RCP<MapType> ghostedSourceMap;
272  {
273  std::vector<GO> indices;
274  sourceGlobalIndexer.getOwnedAndGhostedIndices(indices);
275  ghostedSourceMap = rcp(new MapType(Teuchos::OrdinalTraits<GO>::invalid(),indices,0,comm_));
276  }
277 
278  RCP<GraphType> ghostedGraph = rcp(new GraphType(ghostedTargetMap,ghostedSourceMap,0));
279 
280  // Now insert the non-zero pattern per row
281  std::vector<std::string> elementBlockIds;
282  targetGlobalIndexer_->getElementBlockIds(elementBlockIds);
283  std::vector<std::string>::const_iterator blockItr;
284  for (blockItr=elementBlockIds.begin();blockItr!=elementBlockIds.end();++blockItr) {
285  std::string blockId = *blockItr;
286  const std::vector<LO> & elements = targetGlobalIndexer_->getElementBlock(blockId);
287 
288  std::vector<GO> row_gids;
289  std::vector<GO> col_gids;
290 
291  for(std::size_t elmt=0;elmt<elements.size();elmt++) {
292  targetGlobalIndexer_->getElementGIDs(elements[elmt],row_gids);
293  sourceGlobalIndexer.getElementGIDs(elements[elmt],col_gids);
294  for(std::size_t row=0;row<row_gids.size();row++)
295  ghostedGraph->insertGlobalIndices(row_gids[row],col_gids);
296  }
297  }
298 
299  // Fill complete with range and domain map
300  ghostedGraph->fillComplete(ghostedSourceMap,ghostedTargetMap);
301 
302  // *****************
303  // Build owned graph
304  // *****************
305  RCP<MapType> ownedTargetMap;
306  {
307  std::vector<GO> indices;
308  targetGlobalIndexer_->getOwnedIndices(indices);
309  ownedTargetMap = rcp(new MapType(Teuchos::OrdinalTraits<GO>::invalid(),indices,0,comm_));
310  }
311 
312  RCP<const MapType> ownedSourceMap = inputOwnedSourceMap;
313  if (is_null(ownedSourceMap)) {
314  std::vector<GO> indices;
315  sourceGlobalIndexer.getOwnedIndices(indices);
316  ownedSourceMap = rcp(new MapType(Teuchos::OrdinalTraits<GO>::invalid(),indices,0,comm_));
317  }
318 
319  RCP<GraphType> ownedGraph = rcp(new GraphType(ownedTargetMap,0));
320  RCP<const ExportType> exporter = rcp(new ExportType(ghostedTargetMap,ownedTargetMap));
321  ownedGraph->doExport(*ghostedGraph, *exporter, Tpetra::INSERT);
322  ownedGraph->fillComplete(ownedSourceMap,ownedTargetMap);
323 
324  RCP<MatrixType> ghostedMatrix = rcp(new MatrixType(ghostedGraph));
325  RCP<MatrixType> ownedMatrix = rcp(new MatrixType(ownedGraph));
326  // ghostedMatrix.fillComplete();
327  // ghostedMatrix.resumeFill();
328 
329  ghostedMatrix->setAllToScalar(0.0);
330  ownedMatrix->setAllToScalar(0.0);
331  PHX::Device::fence();
332 
333  // *******************
334  // Fill ghosted matrix
335  // *******************
336  for (const auto& block : elementBlockNames_) {
337 
339  const auto& worksets = worksetContainer_->getWorksets(wd);
340  for (const auto& workset : *worksets) {
341 
342  // Get target basis values: current implementation assumes target basis is HGrad
343  const auto& targetBasisValues = workset.getBasisValues(targetBasisDescriptor_,integrationDescriptor_);
344  const auto& targetWeightedBasis = targetBasisValues.weighted_basis_scalar.get_static_view();
345 
346  // Sources can be any basis
347  const auto& sourceBasisValues = workset.getBasisValues(sourceBasisDescriptor,integrationDescriptor_);
348  Kokkos::View<const double***,PHX::Device> sourceUnweightedScalarBasis;
349  Kokkos::View<const double****,PHX::Device> sourceUnweightedVectorBasis;
350  bool useRankThreeBasis = false; // default to gradient or vector basis
351  if ( (sourceBasisDescriptor.getType() == "HGrad") || (sourceBasisDescriptor.getType() == "Const") ) {
352  if (directionIndex == -1) { // Project dof value
353  sourceUnweightedScalarBasis = sourceBasisValues.basis_scalar.get_static_view();
354  useRankThreeBasis = true;
355  }
356  else { // Project dof gradient of scalar basis
357  sourceUnweightedVectorBasis = sourceBasisValues.grad_basis.get_static_view();
358  }
359  }
360  else { // Project vector value
361  sourceUnweightedVectorBasis = sourceBasisValues.basis_vector.get_static_view();
362  }
363 
364  // Get the element local ids
365  Kokkos::View<LO**,PHX::Device> targetLocalIds("buildRHSMatrix: targetLocalIds", workset.numOwnedCells(),
366  targetGlobalIndexer_->getElementBlockGIDCount(block));
367  Kokkos::View<LO**,PHX::Device> sourceLocalIds("buildRHSMatrix: sourceLocalIds", workset.numOwnedCells(),
368  sourceGlobalIndexer.getElementBlockGIDCount(block));
369  {
370  // Remove the ghosted cell ids or the call to getElementLocalIds will spill array bounds
371  const auto cellLocalIdsNoGhost = Kokkos::subview(workset.cell_local_ids_k,std::make_pair(0,workset.numOwnedCells()));
372  targetGlobalIndexer_->getElementLIDs(cellLocalIdsNoGhost,targetLocalIds);
373  sourceGlobalIndexer.getElementLIDs(cellLocalIdsNoGhost,sourceLocalIds);
374  }
375 
376  // Get the offsets
377  Kokkos::View<LO*,PHX::Device> targetFieldOffsets;
378  {
379  const auto fieldIndex = targetGlobalIndexer_->getFieldNum(targetBasisDescriptor_.getType());
380  const std::vector<LO>& offsets = targetGlobalIndexer_->getGIDFieldOffsets(block,fieldIndex);
381  targetFieldOffsets = Kokkos::View<LO*,PHX::Device>("L2Projection:buildRHS:targetFieldOffsets",offsets.size());
382  const auto hostOffsets = Kokkos::create_mirror_view(targetFieldOffsets);
383  for(size_t i=0; i < offsets.size(); ++i)
384  hostOffsets(i) = offsets[i];
385  Kokkos::deep_copy(targetFieldOffsets,hostOffsets);
386  PHX::Device::fence();
387  }
388 
389  Kokkos::View<LO*,PHX::Device> sourceFieldOffsets;
390  {
391  const auto fieldIndex = sourceGlobalIndexer.getFieldNum(sourceFieldName);
392  const std::vector<LO>& offsets = sourceGlobalIndexer.getGIDFieldOffsets(block,fieldIndex);
393  sourceFieldOffsets = Kokkos::View<LO*,PHX::Device>("L2Projection:buildRHS:sourceFieldOffsets",offsets.size());
394  const auto hostOffsets = Kokkos::create_mirror_view(sourceFieldOffsets);
395  for(size_t i=0; i <offsets.size(); ++i)
396  hostOffsets(i) = offsets[i];
397  Kokkos::deep_copy(sourceFieldOffsets,hostOffsets);
398  PHX::Device::fence();
399  }
400 
401  const auto localMatrix = ghostedMatrix->getLocalMatrix();
402  const int numRows = static_cast<int>(targetWeightedBasis.extent(1));
403  int tmpNumCols = -1;
404  int tmpNumQP = -1;
405  if (useRankThreeBasis) {
406  tmpNumCols = static_cast<int>(sourceUnweightedScalarBasis.extent(1));
407  tmpNumQP = static_cast<int>(sourceUnweightedScalarBasis.extent(2));
408  }
409  else {
410  tmpNumCols = static_cast<int>(sourceUnweightedVectorBasis.extent(1));
411  tmpNumQP = static_cast<int>(sourceUnweightedVectorBasis.extent(2));
412  }
413  const int numCols = tmpNumCols;
414  const int numQP = tmpNumQP;
415 
416  Kokkos::parallel_for(Kokkos::RangePolicy<PHX::Device>(0,workset.numOwnedCells()),KOKKOS_LAMBDA (const int& cell) {
417  LO cLIDs[256];
418  double vals[256];
419  for (int row = 0; row < numRows; ++row) {
420  const int rowOffset = targetFieldOffsets(row);
421  const int rowLID = targetLocalIds(cell,rowOffset);
422  for (int col = 0; col < numCols; ++col)
423  vals[col] = 0.0;
424 
425  for (int col = 0; col < numCols; ++col) {
426  for (int qp = 0; qp < numQP; ++qp) {
427  const int colOffset = sourceFieldOffsets(col);
428  const int colLID = sourceLocalIds(cell,colOffset);
429  cLIDs[col] = colLID;
430  if (useRankThreeBasis)
431  vals[col] += sourceUnweightedScalarBasis(cell,col,qp) * targetWeightedBasis(cell,row,qp);
432  else
433  vals[col] += sourceUnweightedVectorBasis(cell,col,qp,directionIndex) * targetWeightedBasis(cell,row,qp);
434  }
435  }
436  localMatrix.sumIntoValues(rowLID,cLIDs,numCols,vals,true,true);
437  }
438  });
439  }
440 
441  }
442  ghostedMatrix->fillComplete(ghostedSourceMap,ghostedTargetMap);
443  ownedMatrix->resumeFill();
444  ownedMatrix->doExport(*ghostedMatrix,*exporter,Tpetra::ADD);
445  ownedMatrix->fillComplete(ownedSourceMap,ownedTargetMap);
446 
447  return ownedMatrix;
448  }
449 
450 }
451 
452 #endif
is_null
bool is_null(const boost::shared_ptr< T > &p)
Teuchos::MpiComm
Definition: Panzer_L2Projection.hpp:18
panzer::L2Projection::getTargetGlobalIndexer
Teuchos::RCP< panzer::UniqueGlobalIndexer< LO, GO > > getTargetGlobalIndexer() const
Returns the target global indexer. Will be null if setup() has not been called.
Definition: Panzer_L2Projection_impl.hpp:71
Panzer_IntrepidFieldPattern.hpp
Panzer_WorksetContainer.hpp
panzer::BlockedTpetraLinearObjFactory::getGhostedExport
virtual Teuchos::RCP< const ExportType > getGhostedExport(int j) const
get exporter for converting an overalapped object to a "normal" object
Definition: Panzer_BlockedTpetraLinearObjFactory_impl.hpp:931
panzer::IntegrationDescriptor
Definition: Panzer_IntegrationDescriptor.hpp:51
panzer::UniqueGlobalIndexer::getElementGIDs
virtual void getElementGIDs(LocalOrdinalT localElmtId, std::vector< GlobalOrdinalT > &gids, const std::string &blockIdHint="") const =0
Get the global IDs for a particular element. This function overwrites the gids variable.
Panzer_TpetraLinearObjContainer.hpp
panzer::UniqueGlobalIndexer::getOwnedAndGhostedIndices
virtual void getOwnedAndGhostedIndices(std::vector< GlobalOrdinalT > &indices) const =0
Get the set of owned and ghosted indices for this processor.
Panzer_Workset.hpp
panzer::UniqueGlobalIndexer
Definition: Panzer_GatherOrientation_decl.hpp:61
Teuchos::OrdinalTraits
panzer::L2Projection::setup
void setup(const panzer::BasisDescriptor &targetBasis, const panzer::IntegrationDescriptor &integrationDescriptor, const Teuchos::RCP< const Teuchos::MpiComm< int >> &comm, const Teuchos::RCP< const panzer::ConnManager< LO, GO >> &connManager, const std::vector< std::string > &elementBlockNames, const Teuchos::RCP< panzer::WorksetContainer > worksetContainer=Teuchos::null)
Setup base objects for L2 Projections - requires target scalar basis and creates worksets if not supp...
Definition: Panzer_L2Projection_impl.hpp:26
panzer::L2Projection::buildRHSMatrix
Teuchos::RCP< Tpetra::CrsMatrix< double, LO, GO, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device > > > buildRHSMatrix(const panzer::UniqueGlobalIndexer< LO, GO > &sourceDOFManager, const Teuchos::RCP< const Tpetra::Map< LO, GO, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device >>> &ownedSourceMap, const std::string &sourceFieldName, const panzer::BasisDescriptor &sourceBasisDescriptor, const int vectorOrGradientDirectionIndex=-1)
Allocates, fills and returns a rectangular matrix for L2 projection of a scalar field,...
Definition: Panzer_L2Projection_impl.hpp:244
TEUCHOS_ASSERT
#define TEUCHOS_ASSERT(assertion_test)
panzer::BlockedTpetraLinearObjFactory
Definition: Panzer_BlockedTpetraLinearObjFactory.hpp:81
Teuchos::rcp
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
panzer::L2Projection::buildMassMatrix
Teuchos::RCP< Tpetra::CrsMatrix< double, LO, GO, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device > > > buildMassMatrix()
Allocates, fills and returns a mass matrix for L2 projection onto a target basis.
Definition: Panzer_L2Projection_impl.hpp:76
panzer::UniqueGlobalIndexer::getElementLIDs
const Kokkos::View< const LocalOrdinalT *, Kokkos::LayoutRight, PHX::Device > getElementLIDs(LocalOrdinalT localElmtId) const
Definition: Panzer_UniqueGlobalIndexer.hpp:302
panzer::UniqueGlobalIndexer::getGIDFieldOffsets
virtual const std::vector< int > & getGIDFieldOffsets(const std::string &blockId, int fieldNum) const =0
Use the field pattern so that you can find a particular field in the GIDs array.
Teuchos::RCP
panzer::ALL_ELEMENTS
Workset size is set to the total number of local elements in the MPI process.
Definition: Panzer_WorksetDescriptor.hpp:58
panzer::BasisDescriptor
Definition: Panzer_BasisDescriptor.hpp:53
Teuchos_DefaultMpiComm.hpp
panzer::L2Projection::buildInverseLumpedMassMatrix
Teuchos::RCP< Tpetra::MultiVector< double, LO, GO, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device > > > buildInverseLumpedMassMatrix()
Allocates, fills and returns a Tpetra::MultiVector containing the inverse lumped mass matrix values....
Definition: Panzer_L2Projection_impl.hpp:230
panzer::WorksetDescriptor
Definition: Panzer_WorksetDescriptor.hpp:82
panzer::UniqueGlobalIndexer::getElementBlockGIDCount
virtual int getElementBlockGIDCount(const std::string &blockId) const =0
How many GIDs are associate with a particular element block.
panzer::UniqueGlobalIndexer::getFieldNum
virtual int getFieldNum(const std::string &str) const =0
Get the number used for access to this field.
panzer::ConnManager
Definition: Panzer_L2Projection.hpp:25
panzer::UniqueGlobalIndexer::getOwnedIndices
virtual void getOwnedIndices(std::vector< GlobalOrdinalT > &indices) const =0
Get the set of indices owned by this processor.
Panzer_BasisDescriptor.hpp
panzer::Intrepid2FieldPattern
Definition: Panzer_IntrepidFieldPattern.hpp:62
panzer::BlockedTpetraLinearObjFactory::getGhostedTpetraMatrix
Teuchos::RCP< CrsMatrixType > getGhostedTpetraMatrix(int i, int j) const
Definition: Panzer_BlockedTpetraLinearObjFactory_impl.hpp:1064
panzer::BlockedTpetraLinearObjFactory::getTpetraMatrix
Teuchos::RCP< CrsMatrixType > getTpetraMatrix(int i, int j) const
Definition: Panzer_BlockedTpetraLinearObjFactory_impl.hpp:1049
panzer
Definition: Panzer_BasisValues_Evaluator_decl.hpp:54
Panzer_BlockedTpetraLinearObjFactory.hpp
panzer::DOFManager
Definition: Panzer_L2Projection.hpp:26
Panzer_ElementBlockIdToPhysicsIdMap.hpp
panzer::BasisDescriptor::getType
const std::string & getType() const
Get type of basis.
Definition: Panzer_BasisDescriptor.hpp:78
offsets
Kokkos::View< const int *, PHX::Device > offsets
Definition: Panzer_DOFCurl_impl.hpp:180
Panzer_WorksetDescriptor.hpp