Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ScatterDirichletResidual_Epetra_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_SCATTER_DIRICHLET_RESIDUAL_EPETRA_DECL_HPP
12 #define PANZER_EVALUATOR_SCATTER_DIRICHLET_RESIDUAL_EPETRA_DECL_HPP
13 
14 #include "Phalanx_config.hpp"
15 #include "Phalanx_Evaluator_Macros.hpp"
16 #include "Phalanx_MDField.hpp"
17 
19 
20 #include "PanzerDiscFE_config.hpp"
21 #include "Panzer_Dimension.hpp"
22 #include "Panzer_Traits.hpp"
24 
26 
27 class Epetra_Vector;
28 
29 namespace panzer {
30 
31 class EpetraLinearObjContainer;
32 class GlobalIndexer;
33 
42 template<typename EvalT, typename TRAITS,typename LO,typename GO> class ScatterDirichletResidual_Epetra;
43 
44 // **************************************************************
45 // **************************************************************
46 // * Specializations
47 // **************************************************************
48 // **************************************************************
49 
50 
51 // **************************************************************
52 // Residual
53 // **************************************************************
54 template<typename TRAITS,typename LO,typename GO>
56  : public panzer::EvaluatorWithBaseImpl<TRAITS>,
57  public PHX::EvaluatorDerived<panzer::Traits::Residual, TRAITS>,
59 
60 public:
62  const Teuchos::RCP<const panzer::GlobalIndexer> & /* cIndexer=Teuchos::null */)
63  : globalIndexer_(indexer) {}
64 
67  const Teuchos::ParameterList& p);
68 
69  void postRegistrationSetup(typename TRAITS::SetupData d,
71 
72  void preEvaluate(typename TRAITS::PreEvalData d);
73 
74  void evaluateFields(typename TRAITS::EvalData workset);
75 
77  { return Teuchos::rcp(new ScatterDirichletResidual_Epetra<panzer::Traits::Residual,TRAITS,LO,GO>(globalIndexer_,Teuchos::null,pl)); }
78 
79 private:
81 
82  // dummy field so that the evaluator will have something to do
84 
85  // fields that need to be scattered will be put in this vector
86  std::vector< PHX::MDField<const ScalarT,Cell,NODE> > scatterFields_;
87 
88  // maps the local (field,element,basis) triplet to a global ID
89  // for scattering
91  std::vector<int> fieldIds_; // field IDs needing mapping
92 
93  // This maps the scattered field names to the DOF manager field
94  // For instance a Navier-Stokes map might look like
95  // fieldMap_["RESIDUAL_Velocity"] --> "Velocity"
96  // fieldMap_["RESIDUAL_Pressure"] --> "Pressure"
98 
99  std::size_t num_nodes;
100 
101  std::size_t side_subcell_dim_;
102  std::size_t local_side_id_;
103 
105 
107 
108  std::string globalDataKey_; // what global data does this fill?
110 
113 
114  // If set to true, scattering an initial condition
116 
117  // Allows runtime disabling of dirichlet BCs on node-by-node basis
118  std::vector< PHX::MDField<const bool,Cell,NODE> > applyBC_;
119 };
120 
121 // **************************************************************
122 // Tangent
123 // **************************************************************
124 template<typename TRAITS,typename LO,typename GO>
126  : public panzer::EvaluatorWithBaseImpl<TRAITS>,
127  public PHX::EvaluatorDerived<panzer::Traits::Tangent, TRAITS>,
129 
130 public:
132  const Teuchos::RCP<const panzer::GlobalIndexer> & /* cIndexer=Teuchos::null */)
133  : globalIndexer_(indexer) {}
134 
137  const Teuchos::ParameterList& p);
138 
139  void postRegistrationSetup(typename TRAITS::SetupData d,
141 
142  void preEvaluate(typename TRAITS::PreEvalData d);
143 
144  void evaluateFields(typename TRAITS::EvalData workset);
145 
147  { return Teuchos::rcp(new ScatterDirichletResidual_Epetra<panzer::Traits::Tangent,TRAITS,LO,GO>(globalIndexer_,Teuchos::null,pl)); }
148 
149 private:
151 
152  // dummy field so that the evaluator will have something to do
154 
155  // fields that need to be scattered will be put in this vector
156  std::vector< PHX::MDField<const ScalarT,Cell,NODE> > scatterFields_;
157 
158  // maps the local (field,element,basis) triplet to a global ID
159  // for scattering
161  std::vector<int> fieldIds_; // field IDs needing mapping
162 
163  // This maps the scattered field names to the DOF manager field
164  // For instance a Navier-Stokes map might look like
165  // fieldMap_["RESIDUAL_Velocity"] --> "Velocity"
166  // fieldMap_["RESIDUAL_Pressure"] --> "Pressure"
168 
169  std::vector<Teuchos::RCP<Epetra_Vector> > dfdp_vectors_;
170 
171  std::size_t num_nodes;
172 
173  std::size_t side_subcell_dim_;
174  std::size_t local_side_id_;
175 
177 
179 
180  std::string globalDataKey_; // what global data does this fill?
182 
185 
186  // If set to true, scattering an initial condition
188 
189  // Allows runtime disabling of dirichlet BCs on node-by-node basis
190  std::vector< PHX::MDField<const bool,Cell,NODE> > applyBC_;
191 };
192 
193 // **************************************************************
194 // Jacobian
195 // **************************************************************
196 template<typename TRAITS,typename LO,typename GO>
198  : public panzer::EvaluatorWithBaseImpl<TRAITS>,
199  public PHX::EvaluatorDerived<panzer::Traits::Jacobian, TRAITS>,
201 
202 public:
204  const Teuchos::RCP<const panzer::GlobalIndexer> & cIndexer=Teuchos::null)
205  : globalIndexer_(indexer), colGlobalIndexer_(cIndexer) {}
206 
209  const Teuchos::ParameterList& p);
210 
211  void preEvaluate(typename TRAITS::PreEvalData d);
212 
213  void postRegistrationSetup(typename TRAITS::SetupData d,
215 
216  void evaluateFields(typename TRAITS::EvalData workset);
217 
219  { return Teuchos::rcp(new ScatterDirichletResidual_Epetra<panzer::Traits::Jacobian,TRAITS,LO,GO>(globalIndexer_,colGlobalIndexer_,pl)); }
220 
221 private:
222 
224 
225  // dummy field so that the evaluator will have something to do
227 
228  // fields that need to be scattered will be put in this vector
229  std::vector< PHX::MDField<const ScalarT,Cell,NODE> > scatterFields_;
230 
231  // maps the local (field,element,basis) triplet to a global ID
232  // for scattering
234  std::vector<int> fieldIds_; // field IDs needing mapping
235 
236  // This maps the scattered field names to the DOF manager field
237  // For instance a Navier-Stokes map might look like
238  // fieldMap_["RESIDUAL_Velocity"] --> "Velocity"
239  // fieldMap_["RESIDUAL_Pressure"] --> "Pressure"
241 
242  std::size_t num_nodes;
243  std::size_t num_eq;
244 
245  std::size_t side_subcell_dim_;
246  std::size_t local_side_id_;
247 
249 
250  std::string globalDataKey_; // what global data does this fill?
252 
254 
257 
258  // Allows runtime disabling of dirichlet BCs on node-by-node basis
259  std::vector< PHX::MDField<const bool,Cell,NODE> > applyBC_;
260 
262 };
263 
264 }
265 
266 // optionally include hessian support
267 #ifdef Panzer_BUILD_HESSIAN_SUPPORT
269 #endif
270 
271 // **************************************************************
272 #endif
bool checkApplyBC_
If set to true, allows runtime disabling of dirichlet BCs on node-by-node basis.
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
bool checkApplyBC_
If set to true, allows runtime disabling of dirichlet BCs on node-by-node basis.
bool checkApplyBC_
If set to true, allows runtime disabling of dirichlet BCs on node-by-node basis.
ScatterDirichletResidual_Epetra(const Teuchos::RCP< const GlobalIndexer > &indexer, const Teuchos::RCP< const panzer::GlobalIndexer > &cIndexer=Teuchos::null)
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
virtual Teuchos::RCP< CloneableEvaluator > clone(const Teuchos::ParameterList &pl) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
ScatterDirichletResidual_Epetra(const Teuchos::RCP< const GlobalIndexer > &indexer, const Teuchos::RCP< const panzer::GlobalIndexer > &)
ScatterDirichletResidual_Epetra(const Teuchos::RCP< const GlobalIndexer > &indexer, const Teuchos::RCP< const panzer::GlobalIndexer > &)
Non-templated empty base class for template managers.
Pushes residual values into the residual vector for a Newton-based solve.