Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_GatherTangent_Tpetra_decl.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_EVALUATOR_GATHER_TANGENT_TPETRA_DECL_HPP
12 #define PANZER_EVALUATOR_GATHER_TANGENT_TPETRA_DECL_HPP
13 
14 #include "Phalanx_config.hpp"
15 #include "Phalanx_Evaluator_Macros.hpp"
16 #include "Phalanx_MDField.hpp"
17 #include "Phalanx_KokkosViewOfViews.hpp"
18 
20 
21 #include "PanzerDiscFE_config.hpp"
22 #include "Panzer_Dimension.hpp"
23 #include "Panzer_Traits.hpp"
26 
27 #include"Panzer_NodeType.hpp"
28 
30 
31 namespace panzer {
32 
33 class GlobalIndexer; //forward declaration
34 
46 template<typename EvalT,typename TRAITS,typename LO,typename GO,typename NodeT=panzer::TpetraNodeType>
48  : public panzer::EvaluatorWithBaseImpl<TRAITS>,
49  public PHX::EvaluatorDerived<EvalT, TRAITS>,
51 public:
52 
54  globalIndexer_(indexer) {}
55 
57  const Teuchos::ParameterList& p);
58 
59  void postRegistrationSetup(typename TRAITS::SetupData d,
61 
62  void preEvaluate(typename TRAITS::PreEvalData d);
63 
64  void evaluateFields(typename TRAITS::EvalData d);
65 
68 
69  // for testing purposes
70  const PHX::FieldTag & getFieldTag(int i) const
71  { TEUCHOS_ASSERT(i < Teuchos::as<int>(gatherFields_.size())); return gatherFields_[i].fieldTag(); }
72 
73 private:
74 
75  // We always use RealType for gathering as we never compute derivatives for this evaluator
77  //typedef typename EvalT::ScalarT ScalarT;
78 
79  // maps the local (field,element,basis) triplet to a global ID
80  // for scattering
82  std::vector<int> fieldIds_; // field IDs needing mapping
83 
84  std::vector< PHX::MDField<ScalarT,Cell,NODE> > gatherFields_;
86 
89  std::string globalDataKey_; // what global data does this fill?
90 
92 
94 };
95 
96 }
97 
98 // **************************************************************
99 #endif
Teuchos::RCP< std::vector< std::string > > indexerNames_
void preEvaluate(typename TRAITS::PreEvalData d)
void evaluateFields(typename TRAITS::EvalData d)
Teuchos::RCP< const TpetraLinearObjContainer< double, LO, GO, NodeT > > tpetraContainer_
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
const PHX::FieldTag & getFieldTag(int i) const
Gathers tangent vectors dx/dp for computing df/dx*dx/dp + df/dp into the nodal fields of the field ma...
Teuchos::RCP< const panzer::GlobalIndexer > globalIndexer_
std::vector< PHX::MDField< ScalarT, Cell, NODE > > gatherFields_
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
PHX::ViewOfViews< 1, PHX::View< ScalarT ** > > gatherFieldsVoV_
GatherTangent_Tpetra(const Teuchos::RCP< const panzer::GlobalIndexer > &indexer)
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &vm)
#define TEUCHOS_ASSERT(assertion_test)
Non-templated empty base class for template managers.