Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_L2Projection.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef PANZER_L2_PROJECTION_HPP
12 #define PANZER_L2_PROJECTION_HPP
13 
14 #include "Teuchos_RCP.hpp"
15 #include "Phalanx_KokkosDeviceTypes.hpp"
16 #include "PanzerCore_config.hpp"
19 #include "Tpetra_Map.hpp" // for KokkosDeviceWrapperNode
20 #include "Tpetra_CrsMatrix_fwd.hpp"
21 #include "Tpetra_MultiVector_fwd.hpp"
22 #include <vector>
23 #include <string>
24 #include <unordered_map>
25 
26 namespace Teuchos {
27  template<typename T> class MpiComm;
28 }
29 
30 namespace panzer {
31 
32  class BasisDescriptor;
33  class IntegrationDescriptor;
34  class ConnManager;
35  class DOFManager;
36  class GlobalIndexer;
37  class WorksetContainer;
38  template<typename T> class BasisValues2;
39 
43  class L2Projection {
44 
49  std::vector<std::string> elementBlockNames_;
54  std::map<std::string,Teuchos::RCP<panzer::BasisValues2<double>>> basisValues_;
55 
56  public:
57 
59  L2Projection() : setupCalled_(false) {}
60 
70  void setup(const panzer::BasisDescriptor& targetBasis,
71  const panzer::IntegrationDescriptor& integrationDescriptor,
72  const Teuchos::RCP<const Teuchos::MpiComm<int>>& comm,
73  const Teuchos::RCP<const panzer::ConnManager>& connManager,
74  const std::vector<std::string>& elementBlockNames,
75  const Teuchos::RCP<panzer::WorksetContainer> worksetContainer = Teuchos::null);
76 
85  void useBasisValues(const std::map<std::string,Teuchos::RCP<panzer::BasisValues2<double>>>& map_eblock_to_bv);
86 
89 
98  buildMassMatrix(bool use_lumping=false,
99  const std::unordered_map<std::string,double>* elementBlockMultipliers = nullptr);
100 
125 
144  buildRHSMatrix(const panzer::GlobalIndexer& sourceDOFManager,
145  const Teuchos::RCP<const Tpetra::Map<panzer::LocalOrdinal,panzer::GlobalOrdinal,Tpetra::KokkosCompat::KokkosDeviceWrapperNode<PHX::Device>>>& ownedSourceMap,
146  const std::string& sourceFieldName,
147  const panzer::BasisDescriptor& sourceBasisDescriptor,
148  const int vectorOrGradientDirectionIndex = -1);
149  };
150 
151 }
152 
153 #endif
std::vector< std::string > elementBlockNames_
panzer::BasisDescriptor targetBasisDescriptor_
Teuchos::RCP< Tpetra::CrsMatrix< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, Tpetra::KokkosCompat::KokkosDeviceWrapperNode< PHX::Device > > > buildRHSMatrix(const panzer::GlobalIndexer &sourceDOFManager, const Teuchos::RCP< const Tpetra::Map< panzer::LocalOrdinal, panzer::GlobalOrdinal, Tpetra::KokkosCompat::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.
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_
void useBasisValues(const std::map< std::string, Teuchos::RCP< panzer::BasisValues2< double >>> &map_eblock_to_bv)
Override using the panzer::WorksetContainer and instead use the registered BasisValues object...
std::map< std::string, Teuchos::RCP< panzer::BasisValues2< double > > > basisValues_
Teuchos::RCP< Tpetra::MultiVector< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, Tpetra::KokkosCompat::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, Tpetra::KokkosCompat::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. ...