Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_DOFDiv_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_DOF_DIV_IMPL_HPP
12#define PANZER_DOF_DIV_IMPL_HPP
13
18
19namespace panzer {
20
21//**********************************************************************
22template<typename ScalarT,typename ArrayT>
24 PHX::MDField<ScalarT,Cell,IP> dof_div;
25 PHX::MDField<const ScalarT,Cell,Point> dof_value;
26 const ArrayT div_basis;
27 const int num_fields;
28 const int num_points;
29 const int fad_size;
31
32public:
33 using scratch_view = Kokkos::View<ScalarT* ,typename PHX::DevLayout<ScalarT>::type,typename PHX::exec_space::scratch_memory_space,Kokkos::MemoryUnmanaged>;
34
35 EvaluateDOFDiv_withSens(PHX::MDField<ScalarT,Cell,IP> & in_dof_div,
36 PHX::MDField<const ScalarT,Cell,Point> & in_dof_value,
37 const ArrayT & in_div_basis,
38 const bool in_use_shared_memory = false)
39 : dof_div(in_dof_div),
40 dof_value(in_dof_value),
41 div_basis(in_div_basis),
42 num_fields(static_cast<int>(div_basis.extent(1))),
43 num_points(static_cast<int>(div_basis.extent(2))),
44 fad_size(static_cast<int>(Kokkos::dimension_scalar(in_dof_div.get_view()))),
45 use_shared_memory(in_use_shared_memory)
46 {}
47
48 KOKKOS_INLINE_FUNCTION
49 void operator()(const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
50 {
51 const int cell = team.league_rank();
52
54 // Copy reused data into fast scratch space
55 scratch_view dof_values;
56 scratch_view point_values;
57 if (Sacado::IsADType<ScalarT>::value) {
58 dof_values = scratch_view(team.team_shmem(),num_fields,fad_size);
59 point_values = scratch_view(team.team_shmem(),num_points,fad_size);
60 }
61 else {
62 dof_values = scratch_view(team.team_shmem(),num_fields);
63 point_values = scratch_view(team.team_shmem(),num_points);
64 }
65
66 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_fields), [&] (const int& dof) {
67 dof_values(dof) = dof_value(cell,dof);
68 });
69
70 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points), [&] (const int& pt) {
71 point_values(pt) = 0.0;
72 });
73
74 team.team_barrier();
75
76 // Perform contraction
77 for (int dof=0; dof<num_fields; ++dof) {
78 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points), [&] (const int& pt) {
79 point_values(pt) += dof_values(dof) * div_basis(cell,dof,pt);
80 });
81 }
82
83 // Copy to main memory
84 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points), [&] (const int& pt) {
85 dof_div(cell,pt) = point_values(pt);
86 });
87 }
88 else {
89 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points), [&] (const int& pt) {
90 // first initialize to the right thing (prevents over writing with 0)
91 // then loop over one less basis function
92 // ScalarT & div = dof_div(cell,pt);
93 dof_div(cell,pt) = dof_value(cell, 0) * div_basis(cell, 0, pt);
94 for (int bf=1; bf<num_fields; bf++)
95 dof_div(cell,pt) += dof_value(cell, bf) * div_basis(cell, bf, pt);
96 });
97 }
98 }
99 size_t team_shmem_size(int /* team_size */ ) const
100 {
101 if (not use_shared_memory)
102 return 0;
103
104 size_t bytes;
105 if (Sacado::IsADType<ScalarT>::value)
106 bytes = scratch_view::shmem_size(num_fields,fad_size) + scratch_view::shmem_size(num_points,fad_size);
107 else
108 bytes = scratch_view::shmem_size(num_fields) + scratch_view::shmem_size(num_points);
109 return bytes;
110 }
111
112};
113
114//**********************************************************************
115// MOST EVALUATION TYPES
116//**********************************************************************
117
118//**********************************************************************
119template<typename EvalT, typename TRAITS>
121DOFDiv(const Teuchos::ParameterList & p) :
122 use_descriptors_(false),
123 dof_value( p.get<std::string>("Name"),
124 p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
125 basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
126{
127 Teuchos::RCP<const PureBasis> basis
128 = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
129
130 // Verify that this basis supports the div operation
131 TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsDiv(),std::logic_error,
132 "DOFDiv: Basis of type \"" << basis->name() << "\" does not support DIV");
133 TEUCHOS_TEST_FOR_EXCEPTION(!basis->requiresOrientations(),std::logic_error,
134 "DOFDiv: Basis of type \"" << basis->name() << "\" in DOF Div should require orientations. So we are throwing.");
135
136 // build dof_div
137 dof_div = PHX::MDField<ScalarT,Cell,IP>(p.get<std::string>("Div Name"),
138 p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar );
139
140 // add to evaluation graph
141 this->addEvaluatedField(dof_div);
142 this->addDependentField(dof_value);
143
144 std::string n = "DOFDiv: " + dof_div.fieldTag().name() + " ("+PHX::print<EvalT>()+")";
145 this->setName(n);
146}
147
148//**********************************************************************
149template<typename EvalT, typename TRAITS>
151DOFDiv(const PHX::FieldTag & input,
152 const PHX::FieldTag & output,
153 const panzer::BasisDescriptor & bd,
155 : use_descriptors_(true)
156 , bd_(bd)
157 , id_(id)
158 , dof_value(input)
159{
160 TEUCHOS_ASSERT(bd.getType()=="HDiv");
161
162 // build dof_div
163 dof_div = output;
164
165 // add to evaluation graph
166 this->addEvaluatedField(dof_div);
167 this->addDependentField(dof_value);
168
169 std::string n = "DOFDiv: " + dof_div.fieldTag().name() + " ("+PHX::print<EvalT>()+")";
170 this->setName(n);
171}
172
173//**********************************************************************
174template<typename EvalT, typename TRAITS>
176postRegistrationSetup(typename TRAITS::SetupData sd,
178{
179 if(not use_descriptors_)
180 basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
181}
182
183//**********************************************************************
184template<typename EvalT, typename TRAITS>
186evaluateFields(typename TRAITS::EvalData workset)
187{
188 if (workset.num_cells == 0)
189 return;
190
191 const panzer::BasisValues2<double> & basisValues = use_descriptors_ ? this->wda(workset).getBasisValues(bd_,id_)
192 : *this->wda(workset).bases[basis_index];
193
195 Array div_basis = use_descriptors_ ? basisValues.getDivVectorBasis(false) : Array(basisValues.div_basis);
196
197 const bool use_shared_memory = panzer::HP::inst().useSharedMemory<ScalarT>();
198 auto policy = panzer::HP::inst().teamPolicy<ScalarT,PHX::exec_space>(workset.num_cells);
199 auto f = EvaluateDOFDiv_withSens<ScalarT,Array>(dof_div,dof_value,div_basis,use_shared_memory);
200 Kokkos::parallel_for(this->getName(),policy,f);
201}
202
203//**********************************************************************
204
205//**********************************************************************
206// JACOBIAN EVALUATION TYPES
207//**********************************************************************
208
209//**********************************************************************
210template<typename TRAITS>
212DOFDiv(const Teuchos::ParameterList & p) :
213 use_descriptors_(false),
214 dof_value( p.get<std::string>("Name"),
215 p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
216 basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
217{
218 Teuchos::RCP<const PureBasis> basis
219 = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
220
221 // do you specialize because you know where the basis functions are and can
222 // skip a large number of AD calculations?
223 if(p.isType<Teuchos::RCP<const std::vector<int> > >("Jacobian Offsets Vector")) {
224 offsets = *p.get<Teuchos::RCP<const std::vector<int> > >("Jacobian Offsets Vector");
225 accelerate_jacobian = true; // short cut for identity matrix
226 }
227 else
228 accelerate_jacobian = false; // don't short cut for identity matrix
229 accelerate_jacobian = false; // don't short cut for identity matrix
230
231 // Verify that this basis supports the div operation
232 TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsDiv(),std::logic_error,
233 "DOFDiv: Basis of type \"" << basis->name() << "\" does not support DIV");
234 TEUCHOS_TEST_FOR_EXCEPTION(!basis->requiresOrientations(),std::logic_error,
235 "DOFDiv: Basis of type \"" << basis->name() << "\" in DOF Div should require orientations. So we are throwing.");
236
237 // build dof_div
238 dof_div = PHX::MDField<ScalarT,Cell,IP>(p.get<std::string>("Div Name"),
239 p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar );
240
241 // add to evaluation graph
242 this->addEvaluatedField(dof_div);
243 this->addDependentField(dof_value);
244
245 std::string n = "DOFDiv: " + dof_div.fieldTag().name() + " ("+PHX::print<panzer::Traits::Jacobian>()+")";
246 this->setName(n);
247}
248
249//**********************************************************************
250template<typename TRAITS>
252DOFDiv(const PHX::FieldTag & input,
253 const PHX::FieldTag & output,
254 const panzer::BasisDescriptor & bd,
256 : use_descriptors_(true)
257 , bd_(bd)
258 , id_(id)
259 , dof_value(input)
260{
261 TEUCHOS_ASSERT(bd.getType()=="HDiv");
262
263 // build dof_div
264 dof_div = output;
265
266 accelerate_jacobian = false; // don't short cut for identity matrix
267
268 // add to evaluation graph
269 this->addEvaluatedField(dof_div);
270 this->addDependentField(dof_value);
271
272 std::string n = "DOFDiv: " + dof_div.fieldTag().name() + " ("+PHX::print<panzer::Traits::Jacobian>()+")";
273 this->setName(n);
274}
275
276//**********************************************************************
277template<typename TRAITS>
279postRegistrationSetup(typename TRAITS::SetupData sd,
281{
282 this->utils.setFieldData(dof_value,fm);
283 this->utils.setFieldData(dof_div,fm);
284
285 if(not use_descriptors_)
286 basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
287}
288
289template<typename TRAITS>
291evaluateFields(typename TRAITS::EvalData workset)
292{
293 if (workset.num_cells == 0)
294 return;
295
296 const panzer::BasisValues2<double> & basisValues = use_descriptors_ ? this->wda(workset).getBasisValues(bd_,id_)
297 : *this->wda(workset).bases[basis_index];
298
300 Array div_basis = use_descriptors_ ? basisValues.getDivVectorBasis(false) : Array(basisValues.div_basis);
301
302 if(!accelerate_jacobian) {
303 const bool use_shared_memory = panzer::HP::inst().useSharedMemory<ScalarT>();
304 auto policy = panzer::HP::inst().teamPolicy<ScalarT,PHX::exec_space>(workset.num_cells);
305 auto f = EvaluateDOFDiv_withSens<ScalarT,Array>(dof_div,dof_value,div_basis,use_shared_memory);
306 Kokkos::parallel_for(this->getName(),policy,f);
307 return;
308 }
309
310 TEUCHOS_ASSERT(false);
311}
312
313}
314
315#endif
PHX::View< const int * > offsets
PHX::MDField< const ScalarT, Cell, Point > dof_value
const std::string & getType() const
Get type of basis.
PHX::MDField< const Scalar, Cell, BASIS, IP > ConstArray_CellBasisIP
ConstArray_CellBasisIP getDivVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the divergence of a vector basis evaluated at mesh points.
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
PHX::MDField< const ScalarT, Cell, Point > dof_value
void evaluateFields(typename TRAITS::EvalData d)
EvalT::ScalarT ScalarT
DOFDiv(const Teuchos::ParameterList &p)
PHX::MDField< ScalarT, Cell, IP > dof_div
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
Kokkos::View< ScalarT *,typename PHX::DevLayout< ScalarT >::type, typename PHX::exec_space::scratch_memory_space, Kokkos::MemoryUnmanaged > scratch_view
PHX::MDField< const ScalarT, Cell, Point > dof_value
KOKKOS_INLINE_FUNCTION void operator()(const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
EvaluateDOFDiv_withSens(PHX::MDField< ScalarT, Cell, IP > &in_dof_div, PHX::MDField< const ScalarT, Cell, Point > &in_dof_value, const ArrayT &in_div_basis, const bool in_use_shared_memory=false)
PHX::MDField< ScalarT, Cell, IP > dof_div
static HP & inst()
Private ctor.
Kokkos::TeamPolicy< TeamPolicyProperties... > teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
std::vector< std::string >::size_type getBasisIndex(std::string basis_name, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
Returns the index in the workset bases for a particular BasisIRLayout name.