Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ResponseEvaluatorFactory_Probe_impl.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_ResponseEvaluatorFactory_Probe_impl_hpp__
12 #define __Panzer_ResponseEvaluatorFactory_Probe_impl_hpp__
13 
14 #include <string>
15 
16 #include "PanzerDiscFE_config.hpp"
17 
19 #include "Panzer_PhysicsBlock.hpp"
20 #include "Panzer_CellExtreme.hpp"
23 
24 namespace panzer {
25 
26 template <typename EvalT,typename LO,typename GO>
28 buildResponseObject(const std::string & responseName) const
29 {
30  Teuchos::RCP<ResponseBase> response = Teuchos::rcp(new Response_Probe<EvalT>(responseName,comm_,linearObjFactory_));
31  response->setRequiresDirichletAdjustment(applyDirichletToDerivative_);
32 
33  return response;
34 }
35 
36 template <typename EvalT,typename LO,typename GO>
38 buildAndRegisterEvaluators(const std::string & responseName,
40  const panzer::PhysicsBlock & physicsBlock,
41  const Teuchos::ParameterList & user_data) const
42 {
43  using Teuchos::RCP;
44  using Teuchos::rcp;
45 
46  // build scatter evaluator
47  {
49  (globalIndexer_!=Teuchos::null) ? Teuchos::rcp(new ProbeScatter<LO,GO>(globalIndexer_)) : Teuchos::null;
50  std::string field = (fieldName_=="" ? responseName : fieldName_);
51 
52  // Get basis and integration rule associated with field
53  std::vector<panzer::StrPureBasisPair> blockFields = physicsBlock.getProvidedDOFs();
55  for (auto&& v : blockFields) {
56  if (v.first == field) {
57  basis = v.second;
58  break;
59  }
60  }
61  RCP<panzer::IntegrationRule> ir = physicsBlock.getIntegrationRules().at(cubatureDegree_);
62 
63  // build useful evaluator
66  field,
67  fieldComponent_,
68  point_,
69  *ir,
70  basis,
71  globalIndexer_,
72  scatterObj));
73 
74  this->template registerEvaluator<EvalT>(fm, eval);
75 
76  // require last field
77  fm.template requireField<EvalT>(*eval->evaluatedFields()[0]);
78  }
79 }
80 
81 template <typename EvalT,typename LO,typename GO>
84 {
85  if(PHX::print<EvalT>()==PHX::print<panzer::Traits::Residual>() ||
86  PHX::print<EvalT>()==PHX::print<panzer::Traits::Tangent>() ||
87  PHX::print<EvalT>()==PHX::print<panzer::Traits::Jacobian>()
88  )
89  return true;
90 
91  return false;
92 }
93 
94 }
95 
96 #endif
Object that contains information on the physics and discretization of a block of elements with the SA...
virtual Teuchos::RCP< ResponseBase > buildResponseObject(const std::string &responseName) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const std::map< int, Teuchos::RCP< panzer::IntegrationRule > > & getIntegrationRules() const
Returns the unique set of point rules, key is the unique panzer::PointRule::name() ...
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
virtual void buildAndRegisterEvaluators(const std::string &responseName, PHX::FieldManager< panzer::Traits > &fm, const panzer::PhysicsBlock &physicsBlock, const Teuchos::ParameterList &user_data) const
const std::vector< StrPureBasisPair > & getProvidedDOFs() const