Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_ResponseScatterEvaluator_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_RESPONSE_SCATTER_EVALUATOR_EXTREMEVALUE_IMPL_HPP
12#define PANZER_RESPONSE_SCATTER_EVALUATOR_EXTREMEVALUE_IMPL_HPP
13
14#include <iostream>
15#include <string>
16
17#include "PanzerDiscFE_config.hpp"
18
19#include "Phalanx_Evaluator_Macros.hpp"
20#include "Phalanx_MDField.hpp"
21#include "Phalanx_DataLayout_MDALayout.hpp"
22
24#include "Panzer_Dimension.hpp"
26
27#include "Thyra_SpmdVectorBase.hpp"
28#include "Teuchos_ArrayRCP.hpp"
29
30namespace panzer {
31
35template<typename EvalT, typename Traits>
37ResponseScatterEvaluator_ExtremeValue(const std::string & name,
38 const CellData & cd,
39 bool useMax,
40 const Teuchos::RCP<ExtremeValueScatterBase> & extremeValueScatter)
41 : responseName_(name)
42 , scatterObj_(extremeValueScatter)
43 , useMax_(useMax)
44{
45 using Teuchos::RCP;
46 using Teuchos::rcp;
47
48 std::string dummyName = ResponseBase::buildLookupName(name) + " dummy target";
49
50 // build dummy target tag
51 RCP<PHX::DataLayout> dl_dummy = rcp(new PHX::MDALayout<panzer::Dummy>(0));
52 scatterHolder_ = rcp(new PHX::Tag<ScalarT>(dummyName,dl_dummy));
53 this->addEvaluatedField(*scatterHolder_);
54
55 // build dendent field
56 RCP<PHX::DataLayout> dl_cell = rcp(new PHX::MDALayout<panzer::Cell>(cd.numCells()));
57 cellExtremeValue_ = PHX::MDField<const ScalarT,panzer::Cell>(name,dl_cell);
58 this->addDependentField(cellExtremeValue_);
59
60 std::string n = "Extreme Value Response Scatter: " + name;
61 this->setName(n);
62}
63
64template<typename EvalT, typename Traits>
66ResponseScatterEvaluator_ExtremeValue(const std::string & integrandName,
67 const std::string & responseName,
68 const CellData & cd,
69 bool useMax,
70 const Teuchos::RCP<ExtremeValueScatterBase> & extremeValueScatter)
71 : responseName_(responseName)
72 , scatterObj_(extremeValueScatter)
73 , useMax_(useMax)
74{
75 using Teuchos::RCP;
76 using Teuchos::rcp;
77
78 std::string dummyName = ResponseBase::buildLookupName(responseName) + " dummy target";
79
80 // build dummy target tag
81 RCP<PHX::DataLayout> dl_dummy = rcp(new PHX::MDALayout<panzer::Dummy>(0));
82 scatterHolder_ = rcp(new PHX::Tag<ScalarT>(dummyName,dl_dummy));
83 this->addEvaluatedField(*scatterHolder_);
84
85 // build dendent field
86 RCP<PHX::DataLayout> dl_cell = rcp(new PHX::MDALayout<panzer::Cell>(cd.numCells()));
87 cellExtremeValue_ = PHX::MDField<const ScalarT,panzer::Cell>(integrandName,dl_cell);
88 this->addDependentField(cellExtremeValue_);
89
90 std::string n = "Extreme Value Response Scatter: " + responseName;
91 this->setName(n);
92}
93
94template<typename EvalT, typename Traits>
97{
98 // extract linear object container
99 responseObj_ = Teuchos::rcp_dynamic_cast<Response_ExtremeValue<EvalT> >(
100 d.gedc->getDataObject(ResponseBase::buildLookupName(responseName_)),true);
101}
102
103
104template<typename EvalT, typename Traits>
107{
108 auto hostField = Kokkos::create_mirror_view(cellExtremeValue_.get_view());
109 Kokkos::deep_copy(hostField, cellExtremeValue_.get_view());
110 for(index_t i=0;i<d.num_cells;i++) {
111 if(useMax_)
112 responseObj_->value = (responseObj_->value < hostField(i)) ? hostField(i) : responseObj_->value;
113 else
114 responseObj_->value = (responseObj_->value > hostField(i)) ? hostField(i) : responseObj_->value;
115 }
116}
117
118template < >
121{
122 using Teuchos::RCP;
123 using Teuchos::rcp_dynamic_cast;
124 using Thyra::SpmdVectorBase;
125
126 TEUCHOS_ASSERT(scatterObj_!=Teuchos::null);
127
128 // grab local data for inputing
129 Teuchos::ArrayRCP<double> local_dgdx;
130 RCP<SpmdVectorBase<double> > dgdx = rcp_dynamic_cast<SpmdVectorBase<double> >(responseObj_->getGhostedVector());
131 dgdx->getNonconstLocalData(ptrFromRef(local_dgdx));
132 TEUCHOS_ASSERT(!local_dgdx.is_null());
133
134 scatterObj_->scatterDerivative(cellExtremeValue_,d,this->wda,local_dgdx);
135}
136
137}
138
139#endif
Data for determining cell topology and dimensionality.
std::size_t numCells() const
static std::string buildLookupName(const std::string &responseName)
ResponseScatterEvaluator_ExtremeValue(const std::string &name, const CellData &cd, bool useMax, const Teuchos::RCP< ExtremeValueScatterBase > &functionalScatter)
A constructor with concrete arguments instead of a parameter list.
int num_cells
DEPRECATED - use: numCells()
Teuchos::RCP< GlobalEvaluationDataContainer > gedc