4 #ifndef PANZER_L2_PROJECTION_HPP
5 #define PANZER_L2_PROJECTION_HPP
8 #include "Phalanx_KokkosDeviceTypes.hpp"
9 #include "PanzerCore_config.hpp"
12 #include "Tpetra_Map.hpp"
13 #include "Tpetra_CrsMatrix_fwd.hpp"
14 #include "Tpetra_MultiVector_fwd.hpp"
17 #include <unordered_map>
25 class BasisDescriptor;
26 class IntegrationDescriptor;
30 class WorksetContainer;
64 const std::vector<std::string>& elementBlockNames,
79 const std::unordered_map<std::string,double>* elementBlockMultipliers =
nullptr);
125 const Teuchos::RCP<
const Tpetra::Map<panzer::LocalOrdinal,panzer::GlobalOrdinal,Kokkos::Compat::KokkosDeviceWrapperNode<PHX::Device>>>& ownedSourceMap,
126 const std::string& sourceFieldName,
128 const int vectorOrGradientDirectionIndex = -1);
std::vector< std::string > elementBlockNames_
panzer::BasisDescriptor targetBasisDescriptor_
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< panzer::WorksetContainer > worksetContainer_
Teuchos::RCP< const Teuchos::MpiComm< int > > comm_
Teuchos::RCP< panzer::GlobalIndexer > getTargetGlobalIndexer() const
Returns the target global indexer. Will be null if setup() has not been called.
panzer::IntegrationDescriptor integrationDescriptor_
L2Projection()
Constructor.
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, one dimension of gradient (for hgrad basis), or one dimension of a vector field onto the target scalar basis. If you wish to project all values of a vector field or all the gradients of a scalar field, then you will need three separate RHS matrices to form the RHS for each independently. The vectors must be independent Tpetra vectors to solve multiple right hand sides with the linear solver.
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< panzer::DOFManager > targetGlobalIndexer_
Teuchos::RCP< const panzer::ConnManager > connManager_
Unified set of tools for building objects for lumped and consistent L2 projects between bases...
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. ...