Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_ResponseScatterEvaluator_Functional.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_HPP
12#define PANZER_RESPONSE_SCATTER_EVALUATOR_FUNCTIONAL_HPP
13
14#include <iostream>
15#include <string>
16
17#include "PanzerDiscFE_config.hpp"
18#include "Panzer_Dimension.hpp"
19#include "Panzer_CellData.hpp"
23
24#include "Phalanx_Evaluator_Macros.hpp"
25#include "Phalanx_MDField.hpp"
26
28
29namespace panzer {
30
32public:
34
35 virtual void scatterDerivative(const PHX::MDField<const panzer::Traits::Jacobian::ScalarT,panzer::Cell> & cellIntegral,
38 const std::vector<Teuchos::ArrayRCP<double> > & dgdx) const = 0;
39
40#ifdef Panzer_BUILD_HESSIAN_SUPPORT
41 virtual void scatterHessian(const PHX::MDField<const panzer::Traits::Hessian::ScalarT,panzer::Cell> & cellIntegral,
44 const std::vector<Teuchos::ArrayRCP<double> > & d2gdx2) const = 0;
45#endif
46};
47
48template <typename LO,typename GO>
50public:
51 FunctionalScatter(const Teuchos::RCP<const panzer::GlobalIndexer> & globalIndexer)
52 {
53 if(globalIndexer!=Teuchos::null)
54 ugis_.push_back(globalIndexer);
55 }
56
57 FunctionalScatter(const std::vector<Teuchos::RCP<const panzer::GlobalIndexer> > & ugis)
58 : ugis_(ugis) {}
59
60 void scatterDerivative(const PHX::MDField<const panzer::Traits::Jacobian::ScalarT,panzer::Cell> & cellIntegral,
63 const std::vector<Teuchos::ArrayRCP<double> > & dgdx) const;
64
65#ifdef Panzer_BUILD_HESSIAN_SUPPORT
66 void scatterHessian(const PHX::MDField<const panzer::Traits::Hessian::ScalarT,panzer::Cell> & cellIntegral,
69 const std::vector<Teuchos::ArrayRCP<double> > & d2gdx2) const;
70#endif
71
72private:
73
74 std::vector<Teuchos::RCP<const panzer::GlobalIndexer> > ugis_;
75};
76
80template<typename EvalT, typename Traits>
82 public PHX::EvaluatorDerived<EvalT, Traits> {
83public:
84
86 ResponseScatterEvaluator_Functional(const std::string & name,const CellData & cd,
87 const Teuchos::RCP<FunctionalScatterBase> & functionalScatter);
88 ResponseScatterEvaluator_Functional(const std::string & integrandName,const std::string & responseName,const CellData & cd,
89 const Teuchos::RCP<FunctionalScatterBase> & functionalScatter);
90
91 void evaluateFields(typename Traits::EvalData d);
92
93 void preEvaluate(typename Traits::PreEvalData d);
94
95private:
96 typedef typename EvalT::ScalarT ScalarT;
97
98 std::string responseName_;
99 Teuchos::RCP<Response_Functional<EvalT> > responseObj_;
100
101 Teuchos::RCP<PHX::FieldTag> scatterHolder_; // dummy target
102 PHX::MDField<const ScalarT,panzer::Cell> cellIntegral_; // holds cell integrals
103 Teuchos::RCP<FunctionalScatterBase> scatterObj_;
104};
105
106template <typename LO,typename GO>
107void FunctionalScatter<LO,GO>::scatterDerivative(const PHX::MDField<const panzer::Traits::Jacobian::ScalarT,panzer::Cell> & cellIntegral,
110 const std::vector<Teuchos::ArrayRCP<double> > & dgdx) const
111{
112 // for convenience pull out some objects from workset
113 std::string blockId = wda(workset).block_id;
114
115 std::vector<int> blockOffsets;
116 computeBlockOffsets(blockId,ugis_,blockOffsets);
117
118 TEUCHOS_ASSERT(dgdx.size()==ugis_.size());
119
120 auto cellIntegral_h = Kokkos::create_mirror_view(cellIntegral.get_view());
121 Kokkos::deep_copy(cellIntegral_h, cellIntegral.get_view());
122
123 const std::vector<std::size_t> & localCellIds = wda(workset).cell_local_ids;
124
125 for(std::size_t b=0;b<ugis_.size();b++) {
126 int start = blockOffsets[b];
127
128 auto LIDs = ugis_[b]->getLIDs();
129 auto LIDs_h = Kokkos::create_mirror_view(LIDs);
130 Kokkos::deep_copy(LIDs_h, LIDs);
131
132 Teuchos::ArrayRCP<double> dgdx_b = dgdx[b];
133
134 // scatter operation for each cell in workset
135 for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
136 std::size_t cellLocalId = localCellIds[worksetCellIndex];
137
138 // loop over basis functions
139 for(std::size_t i=0;i<LIDs_h.extent(1);i++) {
140 dgdx_b[LIDs_h(cellLocalId, i)] += cellIntegral_h(worksetCellIndex).dx(start+i); // its possible functional is independent of solution value!
141 }
142 }
143 }
144}
145
146#ifdef Panzer_BUILD_HESSIAN_SUPPORT
147template <typename LO,typename GO>
148void FunctionalScatter<LO,GO>::scatterHessian(const PHX::MDField<const panzer::Traits::Hessian::ScalarT,panzer::Cell> & cellIntegral,
151 const std::vector<Teuchos::ArrayRCP<double> > & d2gdx2) const
152{
153 PHX::View<const LO*> LIDs;
154
155 // for convenience pull out some objects from workset
156 std::string blockId = wda(workset).block_id;
157
158 std::vector<int> blockOffsets;
159 computeBlockOffsets(blockId,ugis_,blockOffsets);
160
161 TEUCHOS_ASSERT(d2gdx2.size()==ugis_.size());
162
163 // scatter operation for each cell in workset
164 const std::vector<std::size_t> & localCellIds = wda(workset).cell_local_ids;
165 for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
166 std::size_t cellLocalId = localCellIds[worksetCellIndex];
167
168 for(std::size_t b=0;b<ugis_.size();b++) {
169 int start = blockOffsets[b];
170
171 LIDs = ugis_[b]->getElementLIDs(cellLocalId);
172
173 Teuchos::ArrayRCP<double> d2gdx2_b = d2gdx2[b];
174
175 // loop over basis functions
176 for(std::size_t i=0;i<LIDs.size();i++) {
177 d2gdx2_b[LIDs[i]] += cellIntegral(worksetCellIndex).dx(start+i).dx(0); // its possible functional is independent of solution value!
178 }
179 }
180 }
181}
182#endif
183
184}
185
186#endif
Data for determining cell topology and dimensionality.
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
virtual void scatterDerivative(const PHX::MDField< const panzer::Traits::Jacobian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &dgdx) const =0
virtual void scatterHessian(const PHX::MDField< const panzer::Traits::Hessian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &d2gdx2) const =0
std::vector< Teuchos::RCP< const panzer::GlobalIndexer > > ugis_
FunctionalScatter(const std::vector< Teuchos::RCP< const panzer::GlobalIndexer > > &ugis)
FunctionalScatter(const Teuchos::RCP< const panzer::GlobalIndexer > &globalIndexer)
void scatterDerivative(const PHX::MDField< const panzer::Traits::Jacobian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &dgdx) const
void scatterHessian(const PHX::MDField< const panzer::Traits::Hessian::ScalarT, panzer::Cell > &cellIntegral, panzer::Traits::EvalData workset, WorksetDetailsAccessor &wda, const std::vector< Teuchos::ArrayRCP< double > > &d2gdx2) const
void computeBlockOffsets(const std::string &blockId, const std::vector< Teuchos::RCP< GlobalIndexer > > &ugis, std::vector< int > &blockOffsets)