Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_ResponseEvaluatorFactory_ExtremeValue_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_ExtremeValue_impl_hpp__
12#define __Panzer_ResponseEvaluatorFactory_ExtremeValue_impl_hpp__
13
14#include <string>
15
16#include "PanzerDiscFE_config.hpp"
17
23
24namespace panzer {
25
26template <typename EvalT,typename LO,typename GO>
28buildResponseObject(const std::string & responseName) const
29{
30 Teuchos::RCP<ResponseBase> response = Teuchos::rcp(new Response_ExtremeValue<EvalT>(responseName,comm_,useMax_,linearObjFactory_));
31 response->setRequiresDirichletAdjustment(applyDirichletToDerivative_);
32
33 return response;
34}
35
36template <typename EvalT,typename LO,typename GO>
38buildAndRegisterEvaluators(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
47 // build integration evaluator (integrate over element)
48 if(requiresCellExtreme_) {
49 std::string field = (quadPointField_=="" ? responseName : quadPointField_);
50
51 // build integration rule to use in cell integral
52 RCP<IntegrationRule> ir = rcp(new IntegrationRule(cubatureDegree_,physicsBlock.cellData()));
53
54 Teuchos::ParameterList pl;
55 // add prefix_ to help identify
56 pl.set("Extreme Name",prefix_+field);
57 pl.set("Field Name",field);
58 pl.set("IR",ir);
59 pl.set("Use Max",useMax_);
60
61 Teuchos::RCP<PHX::Evaluator<panzer::Traits> > eval
62 = Teuchos::rcp(new CellExtreme<EvalT,panzer::Traits>(pl));
63
64 this->template registerEvaluator<EvalT>(fm, eval);
65 }
66
67
68 // build scatter evaluator
69 {
70 Teuchos::RCP<ExtremeValueScatterBase> scatterObj =
71 (globalIndexer_!=Teuchos::null) ? Teuchos::rcp(new ExtremeValueScatter<LO,GO>(globalIndexer_)) : Teuchos::null;
72 std::string field = (quadPointField_=="" ? responseName : quadPointField_);
73 field = prefix_+field; // add prefix to help identify
74
75 // build useful evaluator
76 Teuchos::RCP<PHX::Evaluator<panzer::Traits> > eval
78 responseName,
79 physicsBlock.cellData(),
80 useMax_,
81 scatterObj));
82
83 this->template registerEvaluator<EvalT>(fm, eval);
84
85 // require last field
86 fm.template requireField<EvalT>(*eval->evaluatedFields()[0]);
87 }
88}
89
90template <typename EvalT,typename LO,typename GO>
92typeSupported() const
93{
94 if( PHX::print<EvalT>()==PHX::print<panzer::Traits::Residual>() ||
95 PHX::print<EvalT>()==PHX::print<panzer::Traits::Tangent>()
96 )
97 return true;
98
99 if(PHX::print<EvalT>()==PHX::print<panzer::Traits::Jacobian>())
100 return false;
101
102 return false;
103}
104
105}
106
107#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