|
Panzer
Version of the Day
|
Go to the documentation of this file.
43 #ifndef PANZER_EVALUATOR_GATHER_SOLUTION_BLOCKEDTPETRA_DECL_HPP
44 #define PANZER_EVALUATOR_GATHER_SOLUTION_BLOCKEDTPETRA_DECL_HPP
46 #include "Phalanx_config.hpp"
47 #include "Phalanx_Evaluator_Macros.hpp"
48 #include "Phalanx_MDField.hpp"
52 #include "PanzerDiscFE_config.hpp"
61 template <
typename ScalarT,
typename LocalOrdinalT,
typename GlobalOrdinalT,
typename NodeT>
64 template <
typename LocalOrdinalT,
typename GlobalOrdinalT>
67 template <
typename LocalOrdinalT,
typename GlobalOrdinalT>
77 template <
typename EvalT,
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT=panzer::TpetraNodeType>
80 public PHX::EvaluatorDerived<panzer::Traits::Residual, TRAITS>,
98 { std::cout <<
"unspecialized version of \"GatherSolution_BlockedTpetra::evaluateFields\" on \""+PHX::typeAsString<EvalT>()+
"\" should not be used!" << std::endl;
112 template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
115 public PHX::EvaluatorDerived<panzer::Traits::Residual, TRAITS>,
122 : globalIndexer_(indexer) {}
130 void preEvaluate(
typename TRAITS::PreEvalData d);
177 std::vector< std::vector< PHX::MDField<const ScalarT,Cell,NODE> > >
tangentFields_;
191 template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
194 public PHX::EvaluatorDerived<panzer::Traits::Tangent, TRAITS>,
201 : gidIndexer_(indexer) {}
209 void preEvaluate(
typename TRAITS::PreEvalData d);
246 std::vector< std::vector< PHX::MDField<const ScalarT,Cell,NODE> > >
tangentFields_;
254 template <
typename TRAITS,
typename S,
typename LO,
typename GO,
typename NodeT>
257 public PHX::EvaluatorDerived<panzer::Traits::Jacobian, TRAITS>,
262 : globalIndexer_(indexer) {}
270 void preEvaluate(
typename TRAITS::PreEvalData d);
325 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
Teuchos::RCP< const BlockedDOFManager< LO, GO > > globalIndexer_
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
std::string globalDataKey_
Kokkos::View< LO **, PHX::Device > worksetLIDs_
Local indices for unknowns.
Tpetra::CrsGraph< LO, GO, NodeT > CrsGraphType
Teuchos::RCP< const BlockedTpetraLinearObjContainer< S, LO, GO, NodeT > > blockedContainer_
std::vector< Teuchos::RCP< const panzer::UniqueGlobalIndexer< LO, GO > > > fieldGlobalIndexers_
Vector of global indexers, one for each field to gather, respectively.
Tpetra::Export< LO, GO, NodeT > ExportType
panzer::Traits::Tangent::ScalarT ScalarT
Teuchos::RCP< const BlockedDOFManager< LO, GO > > globalIndexer_
Maps the local (field,element,basis) triplet to a global ID for.
Tpetra::Export< LO, GO, NodeT > ExportType
GatherSolution_BlockedTpetra(const Teuchos::RCP< const BlockedDOFManager< LO, GO > > &indexer)
Tpetra::Map< LO, GO, NodeT > MapType
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
std::string globalDataKey_
TRAITS::RealType RealType
Tpetra::Export< LO, GO, NodeT > ExportType
#define TEUCHOS_ASSERT(assertion_test)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
BlockedTpetraLinearObjContainer< S, LO, GO, NodeT > ContainerType
bool disableSensitivities_
Tpetra::CrsMatrix< S, LO, GO, NodeT > CrsMatrixType
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &vm)
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
panzer::Traits::Residual::ScalarT ScalarT
Non-templated empty base class for template managers.
Tpetra::Map< LO, GO, NodeT > MapType
Tpetra::CrsMatrix< S, LO, GO, NodeT > CrsMatrixType
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
GatherSolution_BlockedTpetra(const Teuchos::RCP< const BlockedDOFManager< LO, GO > > &indexer)
BlockedTpetraLinearObjContainer< S, LO, GO, NodeT > ContainerType
std::vector< int > productVectorBlockIndex_
std::vector< std::string > indexerNames_
std::vector< Kokkos::View< int *, PHX::Device > > fieldOffsets_
Offset into the cell lids for each field.
Gathers solution values from the Newton solution vector into the nodal fields of the field manager.
Kokkos::View< LO **, PHX::Device > worksetLIDs_
Local indices for unknowns.
Teuchos::RCP< const BlockedDOFManager< LO, GO > > gidIndexer_
Tpetra::Import< LO, GO, NodeT > ImportType
Tpetra::CrsMatrix< S, LO, GO, NodeT > CrsMatrixType
bool useTimeDerivativeSolutionVector_
Tpetra::Import< LO, GO, NodeT > ImportType
std::vector< int > fieldIds_
Field IDs in the local product vector block (not global field id)
bool useTimeDerivativeSolutionVector_
std::vector< std::vector< PHX::MDField< const ScalarT, Cell, NODE > > > tangentFields_
Teuchos::RCP< const BlockedTpetraLinearObjContainer< S, LO, GO, NodeT > > blockedContainer_
std::string globalDataKey_
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
std::vector< std::string > indexerNames_
void evaluateFields(typename TRAITS::EvalData d)
GatherSolution_BlockedTpetra(const Teuchos::RCP< const BlockedDOFManager< LO, GO > > &indexer)
Teuchos::RCP< const BlockedTpetraLinearObjContainer< S, LO, GO, NodeT > > blockedContainer_
panzer::Traits::Jacobian EvalT
Tpetra::CrsGraph< LO, GO, NodeT > CrsGraphType
panzer::Traits::Tangent EvalT
Tpetra::Vector< S, LO, GO, NodeT > VectorType
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
Tpetra::Vector< S, LO, GO, NodeT > VectorType
Tpetra::Vector< S, LO, GO, NodeT > VectorType
Tpetra::CrsGraph< LO, GO, NodeT > CrsGraphType
std::vector< std::string > indexerNames_
panzer::Traits::Jacobian::ScalarT ScalarT
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
std::vector< int > fieldIds_
std::vector< std::vector< PHX::MDField< const ScalarT, Cell, NODE > > > tangentFields_
std::vector< int > productVectorBlockIndex_
panzer::Traits::Residual EvalT
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
Tpetra::Map< LO, GO, NodeT > MapType
Kokkos::View< LO *, PHX::Device > blockOffsets_
The offset values of the blocked DOFs per element. Size of number of blocks in the product vector + 1...
std::vector< int > fieldIds_
Tpetra::Import< LO, GO, NodeT > ImportType
GatherSolution_BlockedTpetra(const Teuchos::RCP< const BlockedDOFManager< LO, GO > > &indexer)
BlockedTpetraLinearObjContainer< S, LO, GO, NodeT > ContainerType
bool useTimeDerivativeSolutionVector_
std::vector< Kokkos::View< int *, PHX::Device > > fieldOffsets_
Offset into the cell lids for each field. Size of number of fields to scatter.