Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_ResponseScatterEvaluator_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_RESPONSE_SCATTER_EVALUATOR_FUNCTIONAL_IMPL_HPP
12#define PANZER_RESPONSE_SCATTER_EVALUATOR_FUNCTIONAL_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_DefaultProductVector.hpp"
28#include "Thyra_SpmdVectorBase.hpp"
29#include "Thyra_ProductVectorBase.hpp"
30
31#include "Teuchos_ArrayRCP.hpp"
32
33namespace panzer {
34
38template<typename EvalT, typename Traits>
40ResponseScatterEvaluator_Functional(const std::string & name,
41 const CellData & cd,
42 const Teuchos::RCP<FunctionalScatterBase> & functionalScatter)
43 : responseName_(name)
44 , scatterObj_(functionalScatter)
45{
46 using Teuchos::RCP;
47 using Teuchos::rcp;
48
49 std::string dummyName = ResponseBase::buildLookupName(name) + " dummy target";
50
51 // build dummy target tag
52 RCP<PHX::DataLayout> dl_dummy = rcp(new PHX::MDALayout<panzer::Dummy>(0));
53 scatterHolder_ = rcp(new PHX::Tag<ScalarT>(dummyName,dl_dummy));
54 this->addEvaluatedField(*scatterHolder_);
55
56 // build dendent field
57 RCP<PHX::DataLayout> dl_cell = rcp(new PHX::MDALayout<panzer::Cell>(cd.numCells()));
58 cellIntegral_ = PHX::MDField<const ScalarT,panzer::Cell>(name,dl_cell);
59 this->addDependentField(cellIntegral_);
60
61 std::string n = "Functional Response Scatter: " + name;
62 this->setName(n);
63}
64
65template<typename EvalT, typename Traits>
67ResponseScatterEvaluator_Functional(const std::string & integrandName,
68 const std::string & responseName,
69 const CellData & cd,
70 const Teuchos::RCP<FunctionalScatterBase> & functionalScatter)
71 : responseName_(responseName)
72 , scatterObj_(functionalScatter)
73{
74 using Teuchos::RCP;
75 using Teuchos::rcp;
76
77 std::string dummyName = ResponseBase::buildLookupName(responseName) + " dummy target";
78
79 // build dummy target tag
80 RCP<PHX::DataLayout> dl_dummy = rcp(new PHX::MDALayout<panzer::Dummy>(0));
81 scatterHolder_ = rcp(new PHX::Tag<ScalarT>(dummyName,dl_dummy));
82 this->addEvaluatedField(*scatterHolder_);
83
84 // build dendent field
85 RCP<PHX::DataLayout> dl_cell = rcp(new PHX::MDALayout<panzer::Cell>(cd.numCells()));
86 cellIntegral_ = PHX::MDField<const ScalarT,panzer::Cell>(integrandName,dl_cell);
87 this->addDependentField(cellIntegral_);
88
89 std::string n = "Functional Response Scatter: " + responseName;
90 this->setName(n);
91}
92
93template<typename EvalT, typename Traits>
96{
97 // extract linear object container
98 responseObj_ = Teuchos::rcp_dynamic_cast<Response_Functional<EvalT> >(
99 d.gedc->getDataObject(ResponseBase::buildLookupName(responseName_)),true);
100}
101
102
103template<typename EvalT, typename Traits>
106{
107 auto cellIntegral_h = Kokkos::create_mirror_view ( cellIntegral_.get_static_view());
108 Kokkos::deep_copy(cellIntegral_h, cellIntegral_.get_static_view());
109 for(index_t i=0;i<d.num_cells;i++) {
110 responseObj_->value += cellIntegral_h(i);
111 }
112}
113
114template < >
117{
118 using Teuchos::RCP;
119 using Teuchos::rcp_dynamic_cast;
120 using Thyra::SpmdVectorBase;
122
123 TEUCHOS_ASSERT(scatterObj_!=Teuchos::null);
124 TEUCHOS_ASSERT(responseObj_->getGhostedVector()!=Teuchos::null);
125
126 RCP<ProductVectorBase<double> > prod_dgdx = Thyra::castOrCreateNonconstProductVectorBase(responseObj_->getGhostedVector());
127
128 std::vector<Teuchos::ArrayRCP<double> > local_dgdxs;
129 for(int b=0;b<prod_dgdx->productSpace()->numBlocks();b++) {
130 // grab local data for inputing
131 Teuchos::ArrayRCP<double> local_dgdx;
132 RCP<SpmdVectorBase<double> > dgdx = rcp_dynamic_cast<SpmdVectorBase<double> >(prod_dgdx->getNonconstVectorBlock(b));
133 dgdx->getNonconstLocalData(ptrFromRef(local_dgdx));
134
135 TEUCHOS_ASSERT(!local_dgdx.is_null());
136
137 local_dgdxs.push_back(local_dgdx);
138 }
139
140 scatterObj_->scatterDerivative(cellIntegral_,d,this->wda,local_dgdxs);
141}
142
143#ifdef Panzer_BUILD_HESSIAN_SUPPORT
144template < >
147{
148 using Teuchos::RCP;
149 using Teuchos::rcp_dynamic_cast;
150 using Thyra::SpmdVectorBase;
152
153 TEUCHOS_ASSERT(scatterObj_!=Teuchos::null);
154 TEUCHOS_ASSERT(responseObj_->getGhostedVector()!=Teuchos::null);
155
156 RCP<ProductVectorBase<double> > prod_dgdx = Thyra::castOrCreateNonconstProductVectorBase(responseObj_->getGhostedVector());
157
158 std::vector<Teuchos::ArrayRCP<double> > local_dgdxs;
159 for(int b=0;b<prod_dgdx->productSpace()->numBlocks();b++) {
160 // grab local data for inputing
161 Teuchos::ArrayRCP<double> local_dgdx;
162 RCP<SpmdVectorBase<double> > dgdx = rcp_dynamic_cast<SpmdVectorBase<double> >(prod_dgdx->getNonconstVectorBlock(b));
163 dgdx->getNonconstLocalData(ptrFromRef(local_dgdx));
164
165 TEUCHOS_ASSERT(!local_dgdx.is_null());
166
167 local_dgdxs.push_back(local_dgdx);
168 }
169
170 // TEUCHOS_ASSERT(false);
171 scatterObj_->scatterHessian(cellIntegral_,d,this->wda,local_dgdxs);
172}
173#endif
174
175}
176
177#endif
Data for determining cell topology and dimensionality.
std::size_t numCells() const
static std::string buildLookupName(const std::string &responseName)
ResponseScatterEvaluator_Functional(const std::string &name, const CellData &cd, const Teuchos::RCP< FunctionalScatterBase > &functionalScatter)
A constructor with concrete arguments instead of a parameter list.
int num_cells
DEPRECATED - use: numCells()
Teuchos::RCP< GlobalEvaluationDataContainer > gedc