4 #ifndef PANZER_L2_PROJECTION_IMPL_HPP
5 #define PANZER_L2_PROJECTION_IMPL_HPP
7 #include "Teuchos_DefaultMpiComm.hpp"
8 #include "Tpetra_CrsGraph.hpp"
9 #include "Tpetra_MultiVector.hpp"
10 #include "Tpetra_CrsMatrix.hpp"
14 #include "Panzer_TpetraLinearObjFactory.hpp"
29 const Teuchos::RCP<const panzer::ConnManager>& connManager,
30 const std::vector<std::string>& elementBlockNames,
31 const Teuchos::RCP<panzer::WorksetContainer> worksetContainer)
33 targetBasisDescriptor_ = targetBasis;
34 integrationDescriptor_ = integrationDescriptor;
36 connManager_ = connManager;
37 elementBlockNames_ = elementBlockNames;
38 worksetContainer_ = worksetContainer;
42 targetGlobalIndexer_ =
43 Teuchos::rcp(
new panzer::DOFManager(Teuchos::rcp_const_cast<panzer::ConnManager>(connManager),*(comm->getRawMpiComm())));
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];
56 auto intrepidBasis = panzer::createIntrepid2Basis<PHX::Device,double,double>(targetBasisDescriptor_.getType(),
57 targetBasisDescriptor_.getOrder(),
61 targetGlobalIndexer_->addField(eBlock,targetBasisDescriptor_.getType(),fieldPattern);
64 targetGlobalIndexer_->buildGlobalUnknowns();
69 Teuchos::RCP<panzer::GlobalIndexer>
71 {
return targetGlobalIndexer_;}
73 Teuchos::RCP<Tpetra::CrsMatrix<double,panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>>
75 const std::unordered_map<std::string,double>* elementBlockMultipliers)
77 PANZER_FUNC_TIME_MONITOR_DIFF(
"L2Projection::Build Mass Matrix",TopBuildMassMatrix);
78 TEUCHOS_ASSERT(setupCalled_);
80 if (elementBlockMultipliers !=
nullptr) {
81 TEUCHOS_ASSERT(elementBlockMultipliers->size() == elementBlockNames_.size());
85 std::vector<Teuchos::RCP<const panzer::GlobalIndexer>> indexers;
86 indexers.push_back(targetGlobalIndexer_);
92 ownedMatrix->resumeFill();
93 ownedMatrix->setAllToScalar(0.0);
94 ghostedMatrix->resumeFill();
95 ghostedMatrix->setAllToScalar(0.0);
96 PHX::Device().fence();
98 auto M = ghostedMatrix->getLocalMatrix();
99 const int fieldIndex = targetGlobalIndexer_->getFieldNum(targetBasisDescriptor_.getType());
101 const bool is_scalar = targetBasisDescriptor_.getType()==
"HGrad" || targetBasisDescriptor_.getType()==
"Const" || targetBasisDescriptor_.getType()==
"HVol";
105 for (
const auto& block : elementBlockNames_) {
107 double ebMultiplier = 1.0;
108 if (elementBlockMultipliers !=
nullptr)
109 ebMultiplier = elementBlockMultipliers->find(block)->second;
113 const auto worksets = worksetContainer_->getWorksets(wd);
115 for (
const auto& workset : *worksets) {
117 const auto basisValues = workset.getBasisValues(targetBasisDescriptor_,integrationDescriptor_);
119 const auto unweightedBasis = basisValues.basis_scalar;
120 const auto weightedBasis = basisValues.weighted_basis_scalar;
123 const std::vector<panzer::LocalOrdinal>&
offsets = targetGlobalIndexer_->getGIDFieldOffsets(block,fieldIndex);
124 PHX::View<panzer::LocalOrdinal*> kOffsets(
"MassMatrix: Offsets",
offsets.size());
128 PHX::Device().fence();
131 Kokkos::View<panzer::LocalOrdinal**,PHX::Device> localIds(
"MassMatrix: LocalIds", workset.numOwnedCells()+workset.numGhostCells()+workset.numVirtualCells(),
132 targetGlobalIndexer_->getElementBlockGIDCount(block));
135 const auto cellLocalIdsNoGhost = Kokkos::subview(workset.cell_local_ids_k,std::make_pair(0,workset.numOwnedCells()));
137 targetGlobalIndexer_->getElementLIDs(cellLocalIdsNoGhost,localIds);
139 const int numBasisPoints =
static_cast<int>(weightedBasis.extent(1));
141 Kokkos::parallel_for(workset.numOwnedCells(),KOKKOS_LAMBDA (
const int& cell) {
142 double total_mass = 0.0, trace = 0.0;
144 panzer::LocalOrdinal cLIDs[256];
145 const int numIds =
static_cast<int>(localIds.extent(1));
146 for(
int i=0;i<numIds;++i)
147 cLIDs[i] = localIds(cell,i);
149 double vals[256]={0.0};
150 const int numQP =
static_cast<int>(unweightedBasis.extent(2));
152 for (
int row=0; row < numBasisPoints; ++row) {
153 for (
int col=0; col < numIds; ++col) {
154 for (
int qp=0; qp < numQP; ++qp) {
155 auto tmp = unweightedBasis(cell,row,qp) * weightedBasis(cell,col,qp) * ebMultiplier;
163 for (
int row=0; row < numBasisPoints; ++row) {
164 for (
int col=0; col < numBasisPoints; ++col)
167 int offset = kOffsets(row);
168 panzer::LocalOrdinal lid = localIds(cell,offset);
171 for (
int qp=0; qp < numQP; ++qp)
172 vals[col] += unweightedBasis(cell,row,qp) * weightedBasis(cell,col,qp) * ebMultiplier * total_mass / trace;
174 M.sumIntoValues(lid,cLIDs,numIds,vals,
true,
true);
179 Kokkos::parallel_for(workset.numOwnedCells(),KOKKOS_LAMBDA (
const int& cell) {
181 panzer::LocalOrdinal cLIDs[256];
182 const int numIds =
static_cast<int>(localIds.extent(1));
183 for(
int i=0;i<numIds;++i)
184 cLIDs[i] = localIds(cell,i);
187 const int numQP =
static_cast<int>(unweightedBasis.extent(2));
189 for (
int row=0; row < numBasisPoints; ++row) {
190 int offset = kOffsets(row);
191 panzer::LocalOrdinal lid = localIds(cell,offset);
193 for (
int col=0; col < numIds; ++col) {
195 for (
int qp=0; qp < numQP; ++qp)
196 vals[col] += unweightedBasis(cell,row,qp) * weightedBasis(cell,col,qp) * ebMultiplier;
198 M.sumIntoValues(lid,cLIDs,numIds,vals,
true,
true);
207 for (
const auto& block : elementBlockNames_) {
209 double ebMultiplier = 1.0;
210 if (elementBlockMultipliers !=
nullptr)
211 ebMultiplier = elementBlockMultipliers->find(block)->second;
215 const auto& worksets = worksetContainer_->getWorksets(wd);
217 for (
const auto& workset : *worksets) {
219 const auto basisValues = workset.getBasisValues(targetBasisDescriptor_,integrationDescriptor_);
221 const auto unweightedBasis = basisValues.basis_vector;
222 const auto weightedBasis = basisValues.weighted_basis_vector;
225 const std::vector<panzer::LocalOrdinal>&
offsets = targetGlobalIndexer_->getGIDFieldOffsets(block,fieldIndex);
226 PHX::View<panzer::LocalOrdinal*> kOffsets(
"MassMatrix: Offsets",
offsets.size());
230 PHX::Device().fence();
233 Kokkos::View<panzer::LocalOrdinal**,PHX::Device> localIds(
"MassMatrix: LocalIds", workset.numOwnedCells()+workset.numGhostCells()+workset.numVirtualCells(),
234 targetGlobalIndexer_->getElementBlockGIDCount(block));
237 const auto cellLocalIdsNoGhost = Kokkos::subview(workset.cell_local_ids_k,std::make_pair(0,workset.numOwnedCells()));
239 targetGlobalIndexer_->getElementLIDs(cellLocalIdsNoGhost,localIds);
241 const int numBasisPoints =
static_cast<int>(weightedBasis.extent(1));
242 Kokkos::parallel_for(workset.numOwnedCells(),KOKKOS_LAMBDA (
const int& cell) {
244 panzer::LocalOrdinal cLIDs[256];
245 const int numIds =
static_cast<int>(localIds.extent(1));
246 for(
int i=0;i<numIds;++i)
247 cLIDs[i] = localIds(cell,i);
250 const int numQP =
static_cast<int>(unweightedBasis.extent(2));
252 for (
int qp=0; qp < numQP; ++qp) {
253 for (
int row=0; row < numBasisPoints; ++row) {
254 int offset = kOffsets(row);
255 panzer::LocalOrdinal lid = localIds(cell,offset);
257 for (
int col=0; col < numIds; ++col){
259 for(
int dim=0; dim < static_cast<int>(weightedBasis.extent(3)); ++dim)
260 vals[col] += unweightedBasis(cell,row,qp,dim) * weightedBasis(cell,col,qp,dim) * ebMultiplier;
263 M.sumIntoValues(lid,cLIDs,numIds,vals,
true,
true);
271 PHX::exec_space().fence();
274 PANZER_FUNC_TIME_MONITOR_DIFF(
"Exporting of mass matrix",ExportMM);
275 auto map = factory.
getMap(0);
276 ghostedMatrix->fillComplete(map,map);
278 ownedMatrix->doExport(*ghostedMatrix, *exporter, Tpetra::ADD);
279 ownedMatrix->fillComplete();
284 Teuchos::RCP<Tpetra::MultiVector<double,panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>>
287 PANZER_FUNC_TIME_MONITOR_DIFF(
"L2Projection<panzer::LocalOrdinal,panzer::GlobalOrdinal>::buildInverseLumpedMassMatrix",buildInvLMM);
289 const auto massMatrix = this->buildMassMatrix(
true);
290 const auto lumpedMassMatrix = rcp(
new Tpetra::MultiVector<double,panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>(massMatrix->getDomainMap(),1,
true));
291 const auto tmp = rcp(
new Tpetra::MultiVector<double,panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>(massMatrix->getRangeMap(),1,
false));
294 PANZER_FUNC_TIME_MONITOR_DIFF(
"Apply",Apply);
295 massMatrix->apply(*tmp,*lumpedMassMatrix);
298 PANZER_FUNC_TIME_MONITOR_DIFF(
"Reciprocal",Reciprocal);
299 lumpedMassMatrix->reciprocal(*lumpedMassMatrix);
301 return lumpedMassMatrix;
304 Teuchos::RCP<Tpetra::CrsMatrix<double,panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>>
306 const Teuchos::RCP<
const Tpetra::Map<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>>& inputOwnedSourceMap,
307 const std::string& sourceFieldName,
309 const int directionIndex)
316 using MapType = Tpetra::Map<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
317 using GraphType = Tpetra::CrsGraph<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
318 using ExportType = Tpetra::Export<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
319 using MatrixType = Tpetra::CrsMatrix<double,panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
325 RCP<MapType> ghostedTargetMap;
327 std::vector<panzer::GlobalOrdinal> indices;
328 targetGlobalIndexer_->getOwnedAndGhostedIndices(indices);
329 ghostedTargetMap = rcp(
new MapType(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),indices,0,comm_));
332 RCP<MapType> ghostedSourceMap;
334 std::vector<panzer::GlobalOrdinal> indices;
336 ghostedSourceMap = rcp(
new MapType(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),indices,0,comm_));
342 std::vector<size_t> nEntriesPerRow(ghostedTargetMap->getNodeNumElements(),0);
343 std::vector<std::string> elementBlockIds;
344 targetGlobalIndexer_->getElementBlockIds(elementBlockIds);
345 std::vector<std::string>::const_iterator blockItr;
346 for (blockItr=elementBlockIds.begin();blockItr!=elementBlockIds.end();++blockItr) {
347 std::string blockId = *blockItr;
348 const std::vector<panzer::LocalOrdinal> & elements = targetGlobalIndexer_->getElementBlock(blockId);
350 std::vector<panzer::GlobalOrdinal> row_gids;
351 std::vector<panzer::GlobalOrdinal> col_gids;
353 for(std::size_t elmt=0;elmt<elements.size();elmt++) {
354 targetGlobalIndexer_->getElementGIDs(elements[elmt],row_gids);
356 for(std::size_t row=0;row<row_gids.size();row++) {
357 panzer::LocalOrdinal lid =
358 ghostedTargetMap->getLocalElement(row_gids[row]);
359 nEntriesPerRow[lid] += col_gids.size();
364 Teuchos::ArrayView<const size_t> nEntriesPerRowView(nEntriesPerRow);
365 RCP<GraphType> ghostedGraph = rcp(
new GraphType(ghostedTargetMap,ghostedSourceMap,nEntriesPerRowView,Tpetra::StaticProfile));
367 for (blockItr=elementBlockIds.begin();blockItr!=elementBlockIds.end();++blockItr) {
368 std::string blockId = *blockItr;
369 const std::vector<panzer::LocalOrdinal> & elements = targetGlobalIndexer_->getElementBlock(blockId);
371 std::vector<panzer::GlobalOrdinal> row_gids;
372 std::vector<panzer::GlobalOrdinal> col_gids;
374 for(std::size_t elmt=0;elmt<elements.size();elmt++) {
375 targetGlobalIndexer_->getElementGIDs(elements[elmt],row_gids);
377 for(std::size_t row=0;row<row_gids.size();row++)
378 ghostedGraph->insertGlobalIndices(row_gids[row],col_gids);
382 RCP<MapType> ownedTargetMap;
384 std::vector<panzer::GlobalOrdinal> indices;
385 targetGlobalIndexer_->getOwnedIndices(indices);
386 ownedTargetMap = rcp(
new MapType(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),indices,0,comm_));
389 RCP<const MapType> ownedSourceMap = inputOwnedSourceMap;
390 if (is_null(ownedSourceMap)) {
391 std::vector<panzer::GlobalOrdinal> indices;
393 ownedSourceMap = rcp(
new MapType(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),indices,0,comm_));
397 ghostedGraph->fillComplete(ownedSourceMap,ownedTargetMap);
403 RCP<GraphType> ownedGraph = rcp(
new GraphType(ownedTargetMap,0));
404 RCP<const ExportType> exporter = rcp(
new ExportType(ghostedTargetMap,ownedTargetMap));
405 ownedGraph->doExport(*ghostedGraph, *exporter, Tpetra::INSERT);
406 ownedGraph->fillComplete(ownedSourceMap,ownedTargetMap);
408 RCP<MatrixType> ghostedMatrix = rcp(
new MatrixType(ghostedGraph));
409 RCP<MatrixType> ownedMatrix = rcp(
new MatrixType(ownedGraph));
413 ghostedMatrix->setAllToScalar(0.0);
414 ownedMatrix->setAllToScalar(0.0);
415 PHX::Device().fence();
420 for (
const auto& block : elementBlockNames_) {
423 const auto& worksets = worksetContainer_->getWorksets(wd);
424 for (
const auto& workset : *worksets) {
427 const auto& targetBasisValues = workset.getBasisValues(targetBasisDescriptor_,integrationDescriptor_);
428 const auto& targetWeightedBasis = targetBasisValues.weighted_basis_scalar.get_static_view();
431 const auto& sourceBasisValues = workset.getBasisValues(sourceBasisDescriptor,integrationDescriptor_);
432 Kokkos::View<const double***,PHX::Device> sourceUnweightedScalarBasis;
433 Kokkos::View<const double****,PHX::Device> sourceUnweightedVectorBasis;
434 bool useRankThreeBasis =
false;
435 if ( (sourceBasisDescriptor.
getType() ==
"HGrad") || (sourceBasisDescriptor.
getType() ==
"Const") || (sourceBasisDescriptor.
getType() ==
"HVol") ) {
436 if (directionIndex == -1) {
437 sourceUnweightedScalarBasis = sourceBasisValues.basis_scalar.get_static_view();
438 useRankThreeBasis =
true;
441 sourceUnweightedVectorBasis = sourceBasisValues.grad_basis.get_static_view();
445 sourceUnweightedVectorBasis = sourceBasisValues.basis_vector.get_static_view();
449 Kokkos::View<panzer::LocalOrdinal**,PHX::Device> targetLocalIds(
"buildRHSMatrix: targetLocalIds", workset.numOwnedCells(),
450 targetGlobalIndexer_->getElementBlockGIDCount(block));
451 Kokkos::View<panzer::LocalOrdinal**,PHX::Device> sourceLocalIds(
"buildRHSMatrix: sourceLocalIds", workset.numOwnedCells(),
455 const auto cellLocalIdsNoGhost = Kokkos::subview(workset.cell_local_ids_k,std::make_pair(0,workset.numOwnedCells()));
456 targetGlobalIndexer_->getElementLIDs(cellLocalIdsNoGhost,targetLocalIds);
457 sourceGlobalIndexer.
getElementLIDs(cellLocalIdsNoGhost,sourceLocalIds);
461 Kokkos::View<panzer::LocalOrdinal*,PHX::Device> targetFieldOffsets;
463 const auto fieldIndex = targetGlobalIndexer_->getFieldNum(targetBasisDescriptor_.getType());
464 const std::vector<panzer::LocalOrdinal>&
offsets = targetGlobalIndexer_->getGIDFieldOffsets(block,fieldIndex);
465 targetFieldOffsets = Kokkos::View<panzer::LocalOrdinal*,PHX::Device>(
"L2Projection:buildRHS:targetFieldOffsets",
offsets.size());
466 const auto hostOffsets = Kokkos::create_mirror_view(targetFieldOffsets);
467 for(
size_t i=0; i <
offsets.size(); ++i)
469 Kokkos::deep_copy(targetFieldOffsets,hostOffsets);
470 PHX::Device().fence();
473 Kokkos::View<panzer::LocalOrdinal*,PHX::Device> sourceFieldOffsets;
475 const auto fieldIndex = sourceGlobalIndexer.
getFieldNum(sourceFieldName);
477 TEUCHOS_TEST_FOR_EXCEPTION(
offsets.size() == 0, std::runtime_error,
478 "ERROR: panzer::L2Projection::buildRHSMatrix() - The source field, \""
479 << sourceFieldName <<
"\", does not exist in element block \""
481 sourceFieldOffsets = Kokkos::View<panzer::LocalOrdinal*,PHX::Device>(
"L2Projection:buildRHS:sourceFieldOffsets",
offsets.size());
482 const auto hostOffsets = Kokkos::create_mirror_view(sourceFieldOffsets);
483 for(
size_t i=0; i <
offsets.size(); ++i)
485 Kokkos::deep_copy(sourceFieldOffsets,hostOffsets);
486 PHX::Device().fence();
489 const auto localMatrix = ghostedMatrix->getLocalMatrix();
490 const int numRows =
static_cast<int>(targetWeightedBasis.extent(1));
493 if (useRankThreeBasis) {
494 tmpNumCols =
static_cast<int>(sourceUnweightedScalarBasis.extent(1));
495 tmpNumQP =
static_cast<int>(sourceUnweightedScalarBasis.extent(2));
498 tmpNumCols =
static_cast<int>(sourceUnweightedVectorBasis.extent(1));
499 tmpNumQP =
static_cast<int>(sourceUnweightedVectorBasis.extent(2));
501 const int numCols = tmpNumCols;
502 const int numQP = tmpNumQP;
504 Kokkos::parallel_for(Kokkos::RangePolicy<PHX::Device>(0,workset.numOwnedCells()),KOKKOS_LAMBDA (
const int& cell) {
505 panzer::LocalOrdinal cLIDs[256];
507 for (
int row = 0; row < numRows; ++row) {
508 const int rowOffset = targetFieldOffsets(row);
509 const int rowLID = targetLocalIds(cell,rowOffset);
510 for (
int col = 0; col < numCols; ++col)
513 for (
int col = 0; col < numCols; ++col) {
514 for (
int qp = 0; qp < numQP; ++qp) {
515 const int colOffset = sourceFieldOffsets(col);
516 const int colLID = sourceLocalIds(cell,colOffset);
518 if (useRankThreeBasis)
519 vals[col] += sourceUnweightedScalarBasis(cell,col,qp) * targetWeightedBasis(cell,row,qp);
521 vals[col] += sourceUnweightedVectorBasis(cell,col,qp,directionIndex) * targetWeightedBasis(cell,row,qp);
524 localMatrix.sumIntoValues(rowLID,cLIDs,numCols,vals,
true,
true);
530 ghostedMatrix->fillComplete(ownedSourceMap,ownedTargetMap);
531 ownedMatrix->resumeFill();
532 ownedMatrix->doExport(*ghostedMatrix,*exporter,Tpetra::ADD);
533 ownedMatrix->fillComplete(ownedSourceMap,ownedTargetMap);
Kokkos::View< const int *, PHX::Device > offsets
const std::string & getType() const
Get type of basis.
virtual Teuchos::RCP< const MapType > getMap(int i) const
get the map from the matrix
Teuchos::RCP< CrsMatrixType > getTpetraMatrix(int i, int j) const
Teuchos::RCP< CrsMatrixType > getGhostedTpetraMatrix(int i, int j) const
virtual Teuchos::RCP< const ExportType > getGhostedExport(int j) const
get exporter for converting an overalapped object to a "normal" object
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.
const Kokkos::View< const panzer::LocalOrdinal *, Kokkos::LayoutRight, PHX::Device > getElementLIDs(panzer::LocalOrdinal localElmtId) const
virtual void getElementGIDs(panzer::LocalOrdinal localElmtId, std::vector< panzer::GlobalOrdinal > &gids, const std::string &blockIdHint="") const =0
Get the global IDs for a particular element. This function overwrites the gids variable.
virtual void getOwnedAndGhostedIndices(std::vector< panzer::GlobalOrdinal > &indices) const =0
Get the set of owned and ghosted indices for this processor.
virtual int getElementBlockGIDCount(const std::size_t &blockIndex) const =0
How any GIDs are associate with a particular element block.
virtual void getOwnedIndices(std::vector< panzer::GlobalOrdinal > &indices) const =0
Get the set of indices owned by this processor.
virtual int getFieldNum(const std::string &str) const =0
Get the number used for access to this field.
Teuchos::RCP< Tpetra::CrsMatrix< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device > > > buildMassMatrix(bool use_lumping=false, const std::unordered_map< std::string, double > *elementBlockMultipliers=nullptr)
Allocates, fills and returns a mass matrix for L2 projection onto a target basis.
Teuchos::RCP< panzer::GlobalIndexer > getTargetGlobalIndexer() const
Returns the target global indexer. Will be null if setup() has not been called.
void setup(const panzer::BasisDescriptor &targetBasis, const panzer::IntegrationDescriptor &integrationDescriptor, const Teuchos::RCP< const Teuchos::MpiComm< int >> &comm, const Teuchos::RCP< const panzer::ConnManager > &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...
Teuchos::RCP< Tpetra::MultiVector< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device > > > buildInverseLumpedMassMatrix()
Allocates, fills and returns a Tpetra::MultiVector containing the inverse lumped mass matrix values a...
Teuchos::RCP< Tpetra::CrsMatrix< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device > > > buildRHSMatrix(const panzer::GlobalIndexer &sourceDOFManager, const Teuchos::RCP< const Tpetra::Map< panzer::LocalOrdinal, panzer::GlobalOrdinal, 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,...
@ ALL_ELEMENTS
Workset size is set to the total number of local elements in the MPI process.