Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_DOF_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_IMPL_HPP
12#define PANZER_DOF_IMPL_HPP
13
14#include <algorithm>
21
22#include "Intrepid2_FunctionSpaceTools.hpp"
23
24namespace panzer {
25
26//**********************************************************************
27//* DOF evaluator
28//**********************************************************************
29
30//**********************************************************************
31// MOST EVALUATION TYPES
32//**********************************************************************
33
34//**********************************************************************
35template<typename EvalT, typename TRAITS>
37DOF(const Teuchos::ParameterList & p) :
38 use_descriptors_(false),
39 dof_basis( p.get<std::string>("Name"),
40 p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->functional),
41 basis_name(p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->name())
42{
43 Teuchos::RCP<const PureBasis> basis
44 = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
45 is_vector_basis = basis->isVectorBasis();
46
47 // swap between scalar basis value, or vector basis value
48 if(basis->isScalarBasis()) {
49 dof_ip_scalar = PHX::MDField<ScalarT,Cell,Point>(
50 p.get<std::string>("Name"),
51 p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar);
52 this->addEvaluatedField(dof_ip_scalar);
53 }
54 else if(basis->isVectorBasis()) {
55 dof_ip_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(
56 p.get<std::string>("Name"),
57 p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_vector);
58 this->addEvaluatedField(dof_ip_vector);
59 }
60 else
61 { TEUCHOS_ASSERT(false); }
62
63 this->addDependentField(dof_basis);
64
65 std::string n = "DOF: " + dof_basis.fieldTag().name() + " ("+PHX::print<EvalT>()+")";
66 this->setName(n);
67}
68
69//**********************************************************************
70template<typename EvalT, typename TRAITS>
72DOF(const PHX::FieldTag & input,
73 const PHX::FieldTag & output,
74 const panzer::BasisDescriptor & bd,
76 : use_descriptors_(true)
77 , bd_(bd)
78 , id_(id)
79 , dof_basis(input)
80{
81 TEUCHOS_ASSERT(bd.getType()=="HGrad" || bd.getType()=="HCurl" ||
82 bd.getType()=="HDiv" || bd.getType()=="Const")
83
84 is_vector_basis = (bd.getType()=="HCurl" || bd.getType()=="HDiv");
85
86 // swap between scalar basis value, or vector basis value
87 if(not is_vector_basis) {
88 dof_ip_scalar = output;
89 this->addEvaluatedField(dof_ip_scalar);
90 }
91 else {
92 dof_ip_vector = output;
93 this->addEvaluatedField(dof_ip_vector);
94 }
95
96 this->addDependentField(dof_basis);
97
98 std::string n = "DOF: " + dof_basis.fieldTag().name() + " ("+PHX::print<EvalT>()+")";
99 this->setName(n);
100}
101
102//**********************************************************************
103template<typename EvalT, typename TRAITS>
105postRegistrationSetup(typename TRAITS::SetupData sd,
107{
108 this->utils.setFieldData(dof_basis,fm);
109 if(is_vector_basis)
110 this->utils.setFieldData(dof_ip_vector,fm);
111 else
112 this->utils.setFieldData(dof_ip_scalar,fm);
113
114 // descriptors don't access the basis values in the same way
115 if(not use_descriptors_)
116 basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
117}
118
119//**********************************************************************
120template<typename EvalT, typename TRAITS>
122evaluateFields(typename TRAITS::EvalData workset)
123{
124 const panzer::BasisValues2<double> & basisValues = use_descriptors_ ? this->wda(workset).getBasisValues(bd_,id_)
125 : *this->wda(workset).bases[basis_index];
126
127 const auto policy = panzer::HP::inst().teamPolicy<ScalarT,PHX::exec_space>(workset.num_cells);
128 const bool use_shared_memory = panzer::HP::inst().useSharedMemory<ScalarT>();
129
130 if(is_vector_basis) {
132 Array array = use_descriptors_ ? basisValues.getVectorBasisValues(false) : Array(basisValues.basis_vector);
133 const int spaceDim = array.extent(3);
134 if(spaceDim==3) {
135 dof_functors::EvaluateDOFWithSens_Vector<ScalarT,Array,3> functor(dof_basis.get_static_view(),dof_ip_vector.get_static_view(),array,use_shared_memory);
136 Kokkos::parallel_for(this->getName(),policy,functor);
137 }
138 else {
139 dof_functors::EvaluateDOFWithSens_Vector<ScalarT,Array,2> functor(dof_basis.get_static_view(),dof_ip_vector.get_static_view(),array,use_shared_memory);
140 Kokkos::parallel_for(this->getName(),policy,functor);
141 }
142
143 }
144 else {
146 Array interpolation_array = use_descriptors_ ? basisValues.getBasisValues(false) : Array(basisValues.basis_scalar);
147 dof_functors::EvaluateDOFWithSens_Scalar<ScalarT,Array> functor(dof_basis,dof_ip_scalar,interpolation_array);
148 Kokkos::parallel_for(workset.num_cells,functor);
149 }
150}
151
152//**********************************************************************
153
154//**********************************************************************
155// JACOBIAN EVALUATION TYPES
156//**********************************************************************
157
158//**********************************************************************
159template<typename TRAITS>
161DOF(const Teuchos::ParameterList & p) :
162 use_descriptors_(false),
163 dof_basis( p.get<std::string>("Name"),
164 p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->functional),
165 basis_name(p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->name())
166{
167 Teuchos::RCP<const PureBasis> basis
168 = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
169 is_vector_basis = basis->isVectorBasis();
170
171 if(p.isType<Teuchos::RCP<const std::vector<int> > >("Jacobian Offsets Vector")) {
172 const std::vector<int> & offsets = *p.get<Teuchos::RCP<const std::vector<int> > >("Jacobian Offsets Vector");
173
174 // allocate and copy offsets vector to Kokkos array
175 offsets_array = PHX::View<int*>("offsets",offsets.size());
176 auto offsets_array_h = Kokkos::create_mirror_view(offsets_array);
177 for(std::size_t i=0;i<offsets.size();i++)
178 offsets_array_h(i) = offsets[i];
179 Kokkos::deep_copy(offsets_array, offsets_array_h);
180
181 accelerate_jacobian_enabled = true; // short cut for identity matrix
182
183 // get the sensitivities name that is valid for accelerated jacobians
184 sensitivities_name = true;
185 if (p.isType<std::string>("Sensitivities Name"))
186 sensitivities_name = p.get<std::string>("Sensitivities Name");
187 }
188 else
189 accelerate_jacobian_enabled = false; // don't short cut for identity matrix
190
191 // swap between scalar basis value, or vector basis value
192 if(basis->isScalarBasis()) {
193 dof_ip_scalar = PHX::MDField<ScalarT,Cell,Point>(
194 p.get<std::string>("Name"),
195 p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar);
196 this->addEvaluatedField(dof_ip_scalar);
197 }
198 else if(basis->isVectorBasis()) {
199 dof_ip_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(
200 p.get<std::string>("Name"),
201 p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_vector);
202 this->addEvaluatedField(dof_ip_vector);
203 }
204 else
205 { TEUCHOS_ASSERT(false); }
206
207 this->addDependentField(dof_basis);
208
209 std::string n = "DOF: " + dof_basis.fieldTag().name()
210 + ( accelerate_jacobian_enabled ? " accel_jac " : "slow_jac" )
211 + " ("+PHX::print<panzer::Traits::Jacobian>()+")";
212 this->setName(n);
213}
214
215//**********************************************************************
216template<typename TRAITS>
218DOF(const PHX::FieldTag & input,
219 const PHX::FieldTag & output,
220 const panzer::BasisDescriptor & bd,
222 : use_descriptors_(true)
223 , bd_(bd)
224 , id_(id)
225 , dof_basis(input)
226{
227 TEUCHOS_ASSERT(bd.getType()=="HGrad" || bd.getType()=="HCurl" ||
228 bd.getType()=="HDiv" || bd.getType()=="Const")
229
230 accelerate_jacobian_enabled = false; // don't short cut for identity matrix
231
232 is_vector_basis = (bd.getType()=="HCurl" || bd.getType()=="HDiv");
233
234 // swap between scalar basis value, or vector basis value
235 if(not is_vector_basis) {
236 dof_ip_scalar = output;
237 this->addEvaluatedField(dof_ip_scalar);
238 }
239 else {
240 dof_ip_vector = output;
241 this->addEvaluatedField(dof_ip_vector);
242 }
243
244 this->addDependentField(dof_basis);
245
246 std::string n = "DOF: " + dof_basis.fieldTag().name() + " slow_jac(descriptor) ("+PHX::print<typename TRAITS::Jacobian>()+")";
247 this->setName(n);
248}
249
250//**********************************************************************
251template<typename TRAITS>
253postRegistrationSetup(typename TRAITS::SetupData sd,
255{
256 this->utils.setFieldData(dof_basis,fm);
257 if(is_vector_basis)
258 this->utils.setFieldData(dof_ip_vector,fm);
259 else
260 this->utils.setFieldData(dof_ip_scalar,fm);
261
262 // descriptors don't access the basis values in the same way
263 if(not use_descriptors_)
264 basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
265}
266
267// **********************************************************************
268template<typename TRAITS>
269void DOF<typename TRAITS::Jacobian, TRAITS>::
270preEvaluate(typename TRAITS::PreEvalData d)
271{
272 // if sensitivities were requrested for this field enable accelerated
273 // jacobian calculations
274 accelerate_jacobian = false;
275 if(accelerate_jacobian_enabled && d.first_sensitivities_name==sensitivities_name) {
276 accelerate_jacobian = true;
277 }
278}
279
280//**********************************************************************
281template<typename TRAITS>
283evaluateFields(typename TRAITS::EvalData workset)
284{
285 const panzer::BasisValues2<double> & basisValues = use_descriptors_ ? this->wda(workset).getBasisValues(bd_,id_)
286 : *this->wda(workset).bases[basis_index];
287
288 if(is_vector_basis) {
289 if(accelerate_jacobian) {
291 Array array = use_descriptors_ ? basisValues.getVectorBasisValues(false) : Array(basisValues.basis_vector);
292 const int spaceDim = array.extent(3);
293 if(spaceDim==3) {
294 dof_functors::EvaluateDOFFastSens_Vector<ScalarT,Array,3> functor(dof_basis,dof_ip_vector,offsets_array,array);
295 Kokkos::parallel_for(workset.num_cells,functor);
296 }
297 else {
298 dof_functors::EvaluateDOFFastSens_Vector<ScalarT,Array,2> functor(dof_basis,dof_ip_vector,offsets_array,array);
299 Kokkos::parallel_for(workset.num_cells,functor);
300 }
301 }
302 else {
303 const bool use_shared_memory = panzer::HP::inst().useSharedMemory<ScalarT>();
304 const auto policy = panzer::HP::inst().teamPolicy<ScalarT,PHX::exec_space>(workset.num_cells);
306 Array array = use_descriptors_ ? basisValues.getVectorBasisValues(false) : Array(basisValues.basis_vector);
307 const int spaceDim = array.extent(3);
308 if(spaceDim==3) {
309 dof_functors::EvaluateDOFWithSens_Vector<ScalarT,Array,3> functor(dof_basis.get_static_view(),dof_ip_vector.get_static_view(),array,use_shared_memory);
310 Kokkos::parallel_for(this->getName(),policy,functor);
311 }
312 else {
313 dof_functors::EvaluateDOFWithSens_Vector<ScalarT,Array,2> functor(dof_basis.get_static_view(),dof_ip_vector.get_static_view(),array,use_shared_memory);
314 Kokkos::parallel_for(this->getName(),policy,functor);
315 }
316 }
317 }
318 else {
320 Array array = use_descriptors_ ? basisValues.getBasisValues(false) : Array(basisValues.basis_scalar);
321 if(accelerate_jacobian) {
322 dof_functors::EvaluateDOFFastSens_Scalar<ScalarT,Array> functor(dof_basis,dof_ip_scalar,offsets_array,array);
323 Kokkos::parallel_for(workset.num_cells,functor);
324 }
325 else {
326 dof_functors::EvaluateDOFWithSens_Scalar<ScalarT,Array> functor(dof_basis,dof_ip_scalar,array);
327 Kokkos::parallel_for(workset.num_cells,functor);
328 }
329 }
330}
331
332}
333
334#endif
PHX::View< const int * > offsets
const std::string & getType() const
Get type of basis.
PHX::MDField< const Scalar, Cell, BASIS, IP > ConstArray_CellBasisIP
Array_CellBasisIP basis_scalar
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 Scalar, Cell, BASIS, IP, Dim > ConstArray_CellBasisIPDim
Array_CellBasisIPDim basis_vector
ConstArray_CellBasisIPDim getVectorBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at mesh points.
PHX::MDField< const ScalarT, Cell, Point > dof_basis
EvalT::ScalarT ScalarT
PHX::MDField< ScalarT, Cell, Point, Dim > dof_ip_vector
bool is_vector_basis
DOF(const Teuchos::ParameterList &p)
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
void evaluateFields(typename TRAITS::EvalData d)
PHX::MDField< ScalarT, Cell, Point > dof_ip_scalar
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.