Panzer
Version of the Day
|
Unified set of tools for building objects for lumped and consistent L2 projects between bases. More...
#include <Panzer_L2Projection.hpp>
Public Member Functions | |
L2Projection () | |
Constructor. More... | |
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 supplied by user. More... | |
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. More... | |
Teuchos::RCP < panzer::GlobalIndexer > | getTargetGlobalIndexer () const |
Returns the target global indexer. Will be null if setup() has not been called. More... | |
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. More... | |
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 as computed via Hinton 1976. More... | |
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. More... | |
Private Attributes | |
panzer::BasisDescriptor | targetBasisDescriptor_ |
panzer::IntegrationDescriptor | integrationDescriptor_ |
Teuchos::RCP< const Teuchos::MpiComm< int > > | comm_ |
Teuchos::RCP< const panzer::ConnManager > | connManager_ |
std::vector< std::string > | elementBlockNames_ |
Teuchos::RCP < panzer::WorksetContainer > | worksetContainer_ |
bool | setupCalled_ |
Teuchos::RCP< panzer::DOFManager > | targetGlobalIndexer_ |
bool | useUserSuppliedBasisValues_ |
std::map< std::string, Teuchos::RCP < panzer::BasisValues2< double > > > | basisValues_ |
Unified set of tools for building objects for lumped and consistent L2 projects between bases.
Definition at line 43 of file Panzer_L2Projection.hpp.
|
inline |
Constructor.
Definition at line 59 of file Panzer_L2Projection.hpp.
void panzer::L2Projection::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 supplied by user.
[in] | basis | (required) Basis that field values will be projected to. |
[in] | integrationDescriptor | (required) Integration order used for the projection. |
[in] | comm | (required) Teuchos MPI communicator used all processes involved in the project. |
[in] | connManger | (required) Connection manager to describe the mesh. |
[in] | elementBlockNames | (required) Names of element blocks in mesh that are involved in the projection. |
[in] | worksetContainer | (optional) If the user has already allocated worksets for the corresponding mesh/element blocks, we can used those instead of reallocating for projection. |
Definition at line 33 of file Panzer_L2Projection.cpp.
void panzer::L2Projection::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.
[in] | bv | BasisValues object setup in lazy evalaution mode with registered orientations if required. |
The intention of this function is to minimize memory use and redundant evaluations by avoiding workset construction and enforcing lazy evaluation in basis values.
Definition at line 77 of file Panzer_L2Projection.cpp.
Teuchos::RCP< panzer::GlobalIndexer > panzer::L2Projection::getTargetGlobalIndexer | ( | ) | const |
Returns the target global indexer. Will be null if setup() has not been called.
Definition at line 90 of file Panzer_L2Projection.cpp.
Teuchos::RCP< Tpetra::CrsMatrix< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, panzer::TpetraNodeType > > panzer::L2Projection::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.
use_lumping | (optional) If set to true, the returned mass matrix is a lumped diagonal mass matrix following Hinton, et al. 1976. |
elementBlockMultipliers | (optional) If non-null, a multiplier will be used for each element block. The elements should be ordered corresponding to commManger block ordering. |
Definition at line 94 of file Panzer_L2Projection.cpp.
Teuchos::RCP< Tpetra::MultiVector< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, panzer::TpetraNodeType > > panzer::L2Projection::buildInverseLumpedMassMatrix | ( | ) |
Allocates, fills and returns a Tpetra::MultiVector containing the inverse lumped mass matrix values as computed via Hinton 1976.
References:
[1] E. Hinton, T. Rock and O. C. Zienkiewicz, "A Note on Mass Lumping and Related Processes in the Finite Element Method," Earthquake Engineering and Structural Dynamics, 4 (1976), 245-249.
[2] T. J. R. Hughes, "The Finite Element Method: Linear Static and Dynamic Finite Element Analysis," (1987), Section 7.3.2.
NOTE: This is currently sligtly inefficient in that it temporarily allocates the consistent mass matrix. This saves significant code duplication at the temporary runtime cost of allocating the mass matrix. The temporary matrix is deallocated on exit from this function so it should not be an issue.
Definition at line 367 of file Panzer_L2Projection.cpp.
Teuchos::RCP< Tpetra::CrsMatrix< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, panzer::TpetraNodeType > > panzer::L2Projection::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.
[in] | (required) | sourceDOFManager The source dof manger object |
[in] | (required) | ownedSourceMap The Tpetra Map for the owned source vector |
[in] | (required) | sourceFieldName The string name of the source field to project |
[in] | (required) | sourceBasisDescriptor The type of the basis for the source field |
[in] | (optional) | vectorOrGradientDirectionIndex For vector fields, this is the vector index to project (x,y, or z component). For scalar fields, the default value of -1 results in a projection of the scalar field value. If set to 0 or greater, it is assumed that the gradient of the HGrad field is projected and that this value is the dimension index for the particular gradient (x, y, or z component). |
Definition at line 387 of file Panzer_L2Projection.cpp.
|
private |
Definition at line 45 of file Panzer_L2Projection.hpp.
|
private |
Definition at line 46 of file Panzer_L2Projection.hpp.
|
private |
Definition at line 47 of file Panzer_L2Projection.hpp.
|
private |
Definition at line 48 of file Panzer_L2Projection.hpp.
|
private |
Definition at line 49 of file Panzer_L2Projection.hpp.
|
mutableprivate |
Definition at line 50 of file Panzer_L2Projection.hpp.
|
private |
Definition at line 51 of file Panzer_L2Projection.hpp.
|
private |
Definition at line 52 of file Panzer_L2Projection.hpp.
|
private |
Definition at line 53 of file Panzer_L2Projection.hpp.
|
private |
Definition at line 54 of file Panzer_L2Projection.hpp.