Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Integrator_GradBasisCrossVector_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_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP
12#define PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP
13
15//
16// Include Files
17//
19
20// Intrepid2
21#include "Intrepid2_FunctionSpaceTools.hpp"
22
23// Kokkos
24#include "Kokkos_ViewFactory.hpp"
25
26// Panzer
30
31namespace panzer
32{
34 //
35 // Constructor
36 //
38 template<typename EvalT, typename Traits>
42 const std::vector<std::string>& resNames,
43 const std::string& vecName,
44 const panzer::BasisIRLayout& basis,
46 const double& multiplier, /* = 1 */
47 const std::vector<std::string>& fmNames, /* =
48 std::vector<std::string>() */
49 const Teuchos::RCP<PHX::DataLayout>& vecDL /* = Teuchos::null */)
50 :
51 evalStyle_(evalStyle),
52 multiplier_(multiplier),
53 numDims_(resNames.size()),
54 numGradDims_(ir.dl_vector->extent(2)),
55 basisName_(basis.name())
56 {
57 using PHX::View;
58 using panzer::BASIS;
59 using panzer::Cell;
61 using panzer::IP;
62 using PHX::DataLayout;
63 using PHX::Device;
64 using PHX::DevLayout;
65 using PHX::MDField;
66 using std::invalid_argument;
67 using std::logic_error;
68 using std::size_t;
69 using std::string;
70 using Teuchos::RCP;
71
72 // Ensure the input makes sense.
73 TEUCHOS_TEST_FOR_EXCEPTION(numDims_ != 3, invalid_argument, "Error: " \
74 "Integrator_GradBasisCrossVector called with the number of residual " \
75 "names not equal to three.")
76 for (const auto& name : resNames)
77 TEUCHOS_TEST_FOR_EXCEPTION(name == "", invalid_argument, "Error: " \
78 "Integrator_GradBasisCrossVector called with an empty residual name.")
79 TEUCHOS_TEST_FOR_EXCEPTION(vecName == "", invalid_argument, "Error: " \
80 "Integrator_GradBasisCrossVector called with an empty vector name.")
81 RCP<const PureBasis> tmpBasis = basis.getBasis();
82 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsGrad(), logic_error,
83 "Error: Integrator_GradBasisCrossVector: Basis of type \""
84 << tmpBasis->name() << "\" does not support the gradient operator.")
85 RCP<DataLayout> tmpVecDL = ir.dl_vector;
86 if (not vecDL.is_null())
87 {
88 tmpVecDL = vecDL;
89 TEUCHOS_TEST_FOR_EXCEPTION(
90 tmpVecDL->extent(2) < ir.dl_vector->extent(2), logic_error,
91 "Error: Integrator_GradBasisCrossVector: Dimension of space " \
92 "exceeds dimension of Vector Data Layout.")
93 TEUCHOS_TEST_FOR_EXCEPTION(numDims_ !=
94 static_cast<int>(vecDL->extent(2)), logic_error, "Error: " \
95 "Integrator_GradBasisCrossVector: The vector must be the same " \
96 "length as the number of residuals.")
97 } // end if (not vecDL.is_null())
98 TEUCHOS_TEST_FOR_EXCEPTION(numGradDims_ > numDims_, logic_error,
99 "Error: Integrator_GradBasisCrossVector: The vector must have at " \
100 "least as many components as there are dimensions in the mesh.")
101
102 // Create the field for the vector-valued function we're integrating.
103 vector_ = MDField<const ScalarT, Cell, IP, Dim>(vecName, tmpVecDL);
104 this->addDependentField(vector_);
105
106 // Create the fields that we're either contributing to or evaluating
107 // (storing).
108 fields_host_.resize(resNames.size());
109 fields_ = OuterView("Integrator_GradBasisCrossVector::fields_", resNames.size());
110 {
111 int i=0;
112 for (const auto& name : resNames)
113 fields_host_[i++] = MDField<ScalarT, Cell, BASIS>(name, basis.functional);
114 } // end loop over resNames
115
116 for (std::size_t i=0; i< fields_.extent(0); ++i) {
117 const auto& field = fields_host_[i];
119 this->addContributedField(field);
120 else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
121 this->addEvaluatedField(field);
122 }
123
124 // Add the dependent field multipliers, if there are any.
125 int i = 0;
126 fieldMults_.resize(fmNames.size());
127 kokkosFieldMults_ = PHX::View<PHX::UnmanagedView<const ScalarT**>*>("GradBasisCrossVector::KokkosFieldMultipliers", fmNames.size());
128 for (const auto& name : fmNames)
129 {
130 fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
131 this->addDependentField(fieldMults_[i - 1]);
132 } // end loop over the field multipliers
133
134 // Set the name of this object.
135 string n("Integrator_GradBasisCrossVector (");
137 n += "CONTRIBUTES";
138 else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
139 n += "EVALUATES";
140 n += "): {";
141 for (size_t j=0; j < fields_host_.size() - 1; ++j)
142 n += resNames[j] + ", ";
143 n += resNames[resNames.size()-1] + "}";
144 this->setName(n);
145 } // end of Constructor
146
148 //
149 // ParameterList Constructor
150 //
152 template<typename EvalT, typename Traits>
155 const Teuchos::ParameterList& p)
156 :
159 p.get<const std::vector<std::string>>("Residual Names"),
160 p.get<std::string>("Vector Name"),
161 (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
162 (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
163 p.get<double>("Multiplier"),
164 p.isType<Teuchos::RCP<const std::vector<std::string>>>
165 ("Field Multipliers") ?
166 (*p.get<Teuchos::RCP<const std::vector<std::string>>>
167 ("Field Multipliers")) : std::vector<std::string>(),
168 p.isType<Teuchos::RCP<PHX::DataLayout>>("Data Layout Vector") ?
169 p.get<Teuchos::RCP<PHX::DataLayout>>("Data Layout Vector") :
170 Teuchos::null)
171 {
172 using Teuchos::ParameterList;
173 using Teuchos::RCP;
174
175 // Ensure that the input ParameterList didn't contain any bogus entries.
176 RCP<ParameterList> validParams = this->getValidParameters();
177 p.validateParameters(*validParams);
178 } // end of ParameterList Constructor
179
181 //
182 // postRegistrationSetup()
183 //
185 template<typename EvalT, typename Traits>
186 void
189 typename Traits::SetupData sd,
191 {
192 using Kokkos::createDynRankView;
194
195 // Get the PHX::Views of the fields.
196 auto fields_host_mirror_ = Kokkos::create_mirror_view(fields_);
197 for (size_t i=0; i < fields_host_.size(); ++i) {
198 fields_host_mirror_(i) = fields_host_[i].get_static_view();
199 }
200 Kokkos::deep_copy(fields_,fields_host_mirror_);
201
202 // Get the PHX::Views of the field multipliers.
203 auto field_mults_host_mirror_ = Kokkos::create_mirror_view(kokkosFieldMults_);
204 for (size_t i=0; i < fieldMults_.size(); ++i)
205 field_mults_host_mirror_(i) = fieldMults_[i].get_static_view();
206 Kokkos::deep_copy(kokkosFieldMults_,field_mults_host_mirror_);
207
208 // Determine the index in the Workset bases for our particular basis name.
209 basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
210 } // end of postRegistrationSetup()
211
213 //
214 // operator()()
215 //
217 template<typename EvalT, typename Traits>
218 template<int NUM_FIELD_MULT>
219 KOKKOS_INLINE_FUNCTION
220 void
223 const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
224 const size_t& cell) const
225 {
227
228 // Initialize the evaluated fields.
229 const int numBases(fields_[0].extent(1)), numQP(vector_.extent(1));
230 if (evalStyle_ == EvaluatorStyle::EVALUATES)
231 for (int dim(0); dim < numDims_; ++dim)
232 for (int basis(0); basis < numBases; ++basis)
233 fields_[dim](cell, basis) = 0.0;
234
235 // The following if-block is for the sake of optimization depending on the
236 // number of field multipliers.
237 ScalarT tmp[3];
238 const int X(0), Y(1), Z(2);
239 if (NUM_FIELD_MULT == 0)
240 {
241 if (numGradDims_ == 1)
242 {
243 for (int qp(0); qp < numQP; ++qp)
244 {
245 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
246 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
247 for (int basis(0); basis < numBases; ++basis)
248 {
249 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
250 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
251 } // end loop over the bases
252 } // end loop over the quadrature points
253 }
254 else if (numGradDims_ == 2)
255 {
256 for (int qp(0); qp < numQP; ++qp)
257 {
258 tmp[X] = multiplier_ * vector_(cell, qp, X);
259 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
260 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
261 for (int basis(0); basis < numBases; ++basis)
262 {
263 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
264 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
265 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
266 tmp[Y] * basis_(cell, basis, qp, X);
267 } // end loop over the bases
268 } // end loop over the quadrature points
269 }
270 else if (numGradDims_ == 3)
271 {
272 for (int qp(0); qp < numQP; ++qp)
273 {
274 tmp[X] = multiplier_ * vector_(cell, qp, X);
275 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
276 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
277 for (int basis(0); basis < numBases; ++basis)
278 {
279 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
280 tmp[Z] * basis_(cell, basis, qp, Y);
281 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
282 tmp[X] * basis_(cell, basis, qp, Z);
283 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
284 tmp[Y] * basis_(cell, basis, qp, X);
285 } // end loop over the bases
286 } // end loop over the quadrature points
287 } // end if (numGradDims_ == something)
288 }
289 else if (NUM_FIELD_MULT == 1)
290 {
291 if (numGradDims_ == 1)
292 {
293 for (int qp(0); qp < numQP; ++qp)
294 {
295 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
296 vector_(cell, qp, Y);
297 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
298 vector_(cell, qp, Z);
299 for (int basis(0); basis < numBases; ++basis)
300 {
301 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
302 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
303 } // end loop over the bases
304 } // end loop over the quadrature points
305 }
306 else if (numGradDims_ == 2)
307 {
308 for (int qp(0); qp < numQP; ++qp)
309 {
310 tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
311 vector_(cell, qp, X);
312 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
313 vector_(cell, qp, Y);
314 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
315 vector_(cell, qp, Z);
316 for (int basis(0); basis < numBases; ++basis)
317 {
318 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
319 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
320 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
321 tmp[Y] * basis_(cell, basis, qp, X);
322 } // end loop over the bases
323 } // end loop over the quadrature points
324 }
325 else if (numGradDims_ == 3)
326 {
327 for (int qp(0); qp < numQP; ++qp)
328 {
329 tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
330 vector_(cell, qp, X);
331 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
332 vector_(cell, qp, Y);
333 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
334 vector_(cell, qp, Z);
335 for (int basis(0); basis < numBases; ++basis)
336 {
337 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
338 tmp[Z] * basis_(cell, basis, qp, Y);
339 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
340 tmp[X] * basis_(cell, basis, qp, Z);
341 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
342 tmp[Y] * basis_(cell, basis, qp, X);
343 } // end loop over the bases
344 } // end loop over the quadrature points
345 } // end if (numGradDims_ == something)
346 }
347 else
348 {
349 const int numFieldMults(kokkosFieldMults_.extent(0));
350 if (numGradDims_ == 1)
351 {
352 for (int qp(0); qp < numQP; ++qp)
353 {
354 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
355 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
356 for (int fm(0); fm < numFieldMults; ++fm)
357 {
358 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
359 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
360 } // end loop over the field multipliers
361 for (int basis(0); basis < numBases; ++basis)
362 {
363 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
364 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
365 } // end loop over the bases
366 } // end loop over the quadrature points
367 }
368 else if (numGradDims_ == 2)
369 {
370 for (int qp(0); qp < numQP; ++qp)
371 {
372 tmp[X] = multiplier_ * vector_(cell, qp, X);
373 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
374 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
375 for (int fm(0); fm < numFieldMults; ++fm)
376 {
377 tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
378 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
379 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
380 } // end loop over the field multipliers
381 for (int basis(0); basis < numBases; ++basis)
382 {
383 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
384 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
385 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
386 tmp[Y] * basis_(cell, basis, qp, X);
387 } // end loop over the bases
388 } // end loop over the quadrature points
389 }
390 else if (numGradDims_ == 3)
391 {
392 for (int qp(0); qp < numQP; ++qp)
393 {
394 tmp[X] = multiplier_ * vector_(cell, qp, X);
395 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
396 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
397 for (int fm(0); fm < numFieldMults; ++fm)
398 {
399 tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
400 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
401 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
402 } // end loop over the field multipliers
403 for (int basis(0); basis < numBases; ++basis)
404 {
405 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
406 tmp[Z] * basis_(cell, basis, qp, Y);
407 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
408 tmp[X] * basis_(cell, basis, qp, Z);
409 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
410 tmp[Y] * basis_(cell, basis, qp, X);
411 } // end loop over the bases
412 } // end loop over the quadrature points
413 } // end if (numGradDims_ == something)
414 } // end if (NUM_FIELD_MULT == something)
415 } // end of operator()()
416
418 //
419 // evaluateFields()
420 //
422 template<typename EvalT, typename Traits>
423 void
426 typename Traits::EvalData workset)
427 {
428 using Kokkos::parallel_for;
429 using Kokkos::RangePolicy;
430
431 // Grab the basis information.
432 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
433
434 // The following if-block is for the sake of optimization depending on the
435 // number of field multipliers. The parallel_fors will loop over the cells
436 // in the Workset and execute operator()() above.
437 if (fieldMults_.size() == 0)
438 parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
439 else if (fieldMults_.size() == 1)
440 parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
441 else
442 parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
443 } // end of evaluateFields()
444
446 //
447 // getValidParameters()
448 //
450 template<typename EvalT, typename TRAITS>
451 Teuchos::RCP<Teuchos::ParameterList>
454 {
457 using PHX::DataLayout;
458 using std::string;
459 using std::vector;
460 using Teuchos::ParameterList;
461 using Teuchos::RCP;
462 using Teuchos::rcp;
463
464 // Create a ParameterList with all the valid keys we support.
465 RCP<ParameterList> p = rcp(new ParameterList);
466
467 RCP<const vector<string>> resNames;
468 p->set("Residual Names", resNames);
469 p->set<string>("Vector Name", "?");
470 RCP<BasisIRLayout> basis;
471 p->set("Basis", basis);
472 RCP<IntegrationRule> ir;
473 p->set("IR", ir);
474 p->set<double>("Multiplier", 1.0);
475 RCP<const vector<string>> fms;
476 p->set("Field Multipliers", fms);
477 RCP<DataLayout> vecDL;
478 p->set("Data Layout Vector", vecDL);
479
480 return p;
481 } // end of getValidParameters()
482
483} // end of namespace panzer
484
485#endif // PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
double multiplier
The scalar multiplier out in front of the integral ( ).
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
PHX::MDField< const ScalarT, Cell, IP, Dim > vector_
A field representing the vector-valued function we're integrating ( ).
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ,...
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
Integrator_GradBasisCrossVector(const panzer::EvaluatorStyle &evalStyle, const std::vector< std::string > &resNames, const std::string &vecName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >(), const Teuchos::RCP< PHX::DataLayout > &vecDL=Teuchos::null)
Main Constructor.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
int numDims_
The number of dimensions associated with the vector.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > fields_host_
The fields to which we'll contribute, or in which we'll store, the result of computing this integral.
int numGradDims_
The number of dimensions associated with the gradient.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
int num_cells
DEPRECATED - use: numCells()
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.
EvaluatorStyle
An indication of how an Evaluator will behave.
This empty struct allows us to optimize operator()() depending on the number of field multipliers.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_