Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_ResponseEvaluatorFactory_Functional_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_Functional_impl_hpp__
12#define __Panzer_ResponseEvaluatorFactory_Functional_impl_hpp__
13
14#include <string>
15
16#include "PanzerDiscFE_config.hpp"
17
20#include "Panzer_Integrator_Scalar.hpp"
25
26namespace panzer {
27
28template <typename EvalT,typename LO,typename GO>
30buildResponseObject(const std::string & responseName) const
31{
32 Teuchos::RCP<ResponseBase> response = Teuchos::rcp(new Response_Functional<EvalT>(responseName,comm_,linearObjFactory_));
33 response->setRequiresDirichletAdjustment(applyDirichletToDerivative_);
34
35 return response;
36}
37
38template <typename EvalT,typename LO,typename GO>
40buildAndRegisterEvaluators(const std::string & responseName,
42 const panzer::PhysicsBlock & physicsBlock,
43 const Teuchos::ParameterList & /* user_data */) const
44{
45 using Teuchos::RCP;
46 using Teuchos::rcp;
47
48
49 // build integration evaluator (integrate over element)
50 if(requiresCellIntegral_) {
51 std::string field = (quadPointField_=="" ? responseName : quadPointField_);
52
53 // build integration rule to use in cell integral
54 RCP<IntegrationRule> ir = rcp(new IntegrationRule(cubatureDegree_,physicsBlock.cellData()));
55
56 Teuchos::ParameterList pl;
57 pl.set("Integral Name", field);
58 pl.set("Integrand Name",field);
59 pl.set("IR",ir);
60
61 Teuchos::RCP<PHX::Evaluator<panzer::Traits> > eval
62 = Teuchos::rcp(new Integrator_Scalar<EvalT,panzer::Traits>(pl));
63
64 this->template registerEvaluator<EvalT>(fm, eval);
65 }
66
67
68 // build scatter evaluator
69 {
70 Teuchos::RCP<FunctionalScatterBase> scatterObj;
71 if(linearObjFactory_!=Teuchos::null) {
72
73 TEUCHOS_ASSERT(linearObjFactory_->getDomainGlobalIndexer()!=Teuchos::null);
74
75 auto ugi = Teuchos::rcp_dynamic_cast<const GlobalIndexer>(linearObjFactory_->getDomainGlobalIndexer());
76 auto bugi = Teuchos::rcp_dynamic_cast<const BlockedDOFManager>(linearObjFactory_->getDomainGlobalIndexer());
77
78 if(ugi!=Teuchos::null) {
79 std::vector<Teuchos::RCP<const GlobalIndexer> > ugis;
80 ugis.push_back(ugi);
81
82 scatterObj = Teuchos::rcp(new FunctionalScatter<LO,GO>(ugis));
83 }
84 else if(bugi!=Teuchos::null) {
85 scatterObj = Teuchos::rcp(new FunctionalScatter<LO,GO>(nc2c_vector(bugi->getFieldDOFManagers())));
86 }
87 else {
88 TEUCHOS_ASSERT(false); // no real global indexer to use
89 }
90 }
91
92 std::string field = (quadPointField_=="" ? responseName : quadPointField_);
93
94 // build useful evaluator
95 Teuchos::RCP<PHX::Evaluator<panzer::Traits> > eval
96 = Teuchos::rcp(new ResponseScatterEvaluator_Functional<EvalT,panzer::Traits>(field,responseName,physicsBlock.cellData(),scatterObj));
97
98 this->template registerEvaluator<EvalT>(fm, eval);
99
100 // require last field
101 fm.template requireField<EvalT>(*eval->evaluatedFields()[0]);
102 }
103}
104
105template <typename EvalT,typename LO,typename GO>
107typeSupported() const
108{
109 if( PHX::print<EvalT>()==PHX::print<panzer::Traits::Residual>() ||
110 PHX::print<EvalT>()==PHX::print<panzer::Traits::Tangent>()
111 )
112 return true;
113
114 if(PHX::print<EvalT>()==PHX::print<panzer::Traits::Jacobian>())
115 return linearObjFactory_!=Teuchos::null;
116
117#ifdef Panzer_BUILD_HESSIAN_SUPPORT
118 if(PHX::print<EvalT>()==PHX::print<panzer::Traits::Hessian>()) {
119 return linearObjFactory_!=Teuchos::null;
120 }
121#endif
122
123 return false;
124}
125
126}
127
128#endif
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
Object that contains information on the physics and discretization of a block of elements with the SA...
const panzer::CellData & cellData() const
virtual void buildAndRegisterEvaluators(const std::string &responseName, PHX::FieldManager< panzer::Traits > &fm, const panzer::PhysicsBlock &physicsBlock, const Teuchos::ParameterList &user_data) const
virtual Teuchos::RCP< ResponseBase > buildResponseObject(const std::string &responseName) const
std::vector< Teuchos::RCP< const GlobalIndexer > > nc2c_vector(const std::vector< Teuchos::RCP< GlobalIndexer > > &ugis)