43 #ifndef PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP
44 #define PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP
53 #include "Intrepid2_FunctionSpaceTools.hpp"
56 #include "Kokkos_ViewFactory.hpp"
70 template<
typename EvalT,
typename Traits>
74 const std::vector<std::string>& resNames,
75 const std::string& vecName,
79 const std::vector<std::string>& fmNames,
81 const Teuchos::RCP<PHX::DataLayout>& vecDL )
85 numDims_(resNames.size()),
86 numGradDims_(ir.dl_vector->extent(2)),
87 basisName_(basis.name())
94 using PHX::DataLayout;
98 using std::invalid_argument;
99 using std::logic_error;
105 TEUCHOS_TEST_FOR_EXCEPTION(
numDims_ != 3, invalid_argument,
"Error: " \
106 "Integrator_GradBasisCrossVector called with the number of residual " \
107 "names not equal to three.")
108 for (
const auto& name : resNames)
109 TEUCHOS_TEST_FOR_EXCEPTION(name ==
"", invalid_argument,
"Error: " \
110 "Integrator_GradBasisCrossVector called with an empty residual name.")
111 TEUCHOS_TEST_FOR_EXCEPTION(vecName ==
"", invalid_argument,
"Error: " \
112 "Integrator_GradBasisCrossVector called with an empty vector name.")
113 RCP<const PureBasis> tmpBasis = basis.
getBasis();
114 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsGrad(), logic_error,
115 "Error: Integrator_GradBasisCrossVector: Basis of type \""
116 << tmpBasis->name() <<
"\" does not support the gradient operator.")
118 if (not vecDL.is_null())
121 TEUCHOS_TEST_FOR_EXCEPTION(
122 tmpVecDL->extent(2) < ir.
dl_vector->extent(2), logic_error,
123 "Error: Integrator_GradBasisCrossVector: Dimension of space " \
124 "exceeds dimension of Vector Data Layout.")
125 TEUCHOS_TEST_FOR_EXCEPTION(
numDims_ !=
126 static_cast<int>(vecDL->extent(2)), logic_error,
"Error: " \
127 "Integrator_GradBasisCrossVector: The vector must be the same " \
128 "length as the number of residuals.")
131 "Error: Integrator_GradBasisCrossVector: The vector must have at " \
132 "least as many components as there are dimensions in the mesh.")
135 vector_ = MDField<const ScalarT, Cell, IP, Dim>(vecName, tmpVecDL);
136 this->addDependentField(
vector_);
140 fields_ = Kokkos::View<PHX::MDField<ScalarT, Cell, BASIS>*>(
"Integrator_GradBasisCrossVector::fields_", resNames.size());
144 for (
const auto& name : resNames)
147 for (std::size_t i=0; i<
fields_.extent(0); ++i) {
150 this->addContributedField(
field);
152 this->addEvaluatedField(
field);
158 typename DevLayout<ScalarT>::type, Device>*>(
159 "GradBasisCrossVector::KokkosFieldMultipliers", fmNames.size());
160 for (
const auto& name : fmNames)
167 string n(
"Integrator_GradBasisCrossVector (");
173 for (
size_t j=0; j <
fields_.extent(0) - 1; ++j)
174 n +=
fields_[j].fieldTag().name() +
", ";
184 template<
typename EvalT,
typename Traits>
187 const Teuchos::ParameterList& p)
191 p.get<const
std::vector<
std::string>>(
"Residual Names"),
192 p.get<
std::string>(
"Vector Name"),
195 p.get<double>(
"Multiplier"),
197 (
"Field Multipliers") ?
199 (
"Field Multipliers")) :
std::vector<
std::string>(),
200 p.isType<
Teuchos::RCP<
PHX::DataLayout>>(
"Data Layout Vector") ?
201 p.get<
Teuchos::RCP<
PHX::DataLayout>>(
"Data Layout Vector") :
204 using Teuchos::ParameterList;
209 p.validateParameters(*validParams);
217 template<
typename EvalT,
typename Traits>
229 for (
size_t i(0); i < fieldMults_.size(); ++i)
230 kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
242 template<
typename EvalT,
typename Traits>
243 template<
int NUM_FIELD_MULT>
244 KOKKOS_INLINE_FUNCTION
249 const size_t& cell)
const
254 const int numBases(fields_[0].extent(1)), numQP(vector_.extent(1));
256 for (
int dim(0); dim < numDims_; ++dim)
257 for (
int basis(0); basis < numBases; ++basis)
258 fields_[dim](cell, basis) = 0.0;
263 const int X(0), Y(1), Z(2);
264 if (NUM_FIELD_MULT == 0)
266 if (numGradDims_ == 1)
268 for (
int qp(0); qp < numQP; ++qp)
270 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
271 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
272 for (
int basis(0); basis < numBases; ++basis)
274 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
275 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
279 else if (numGradDims_ == 2)
281 for (
int qp(0); qp < numQP; ++qp)
283 tmp[X] = multiplier_ * vector_(cell, qp, X);
284 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
285 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
286 for (
int basis(0); basis < numBases; ++basis)
288 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
289 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
290 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
291 tmp[Y] * basis_(cell, basis, qp, X);
295 else if (numGradDims_ == 3)
297 for (
int qp(0); qp < numQP; ++qp)
299 tmp[X] = multiplier_ * vector_(cell, qp, X);
300 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
301 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
302 for (
int basis(0); basis < numBases; ++basis)
304 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
305 tmp[Z] * basis_(cell, basis, qp, Y);
306 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
307 tmp[X] * basis_(cell, basis, qp, Z);
308 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
309 tmp[Y] * basis_(cell, basis, qp, X);
314 else if (NUM_FIELD_MULT == 1)
316 if (numGradDims_ == 1)
318 for (
int qp(0); qp < numQP; ++qp)
320 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
321 vector_(cell, qp, Y);
322 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
323 vector_(cell, qp, Z);
324 for (
int basis(0); basis < numBases; ++basis)
326 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
327 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
331 else if (numGradDims_ == 2)
333 for (
int qp(0); qp < numQP; ++qp)
335 tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
336 vector_(cell, qp, X);
337 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
338 vector_(cell, qp, Y);
339 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
340 vector_(cell, qp, Z);
341 for (
int basis(0); basis < numBases; ++basis)
343 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
344 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
345 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
346 tmp[Y] * basis_(cell, basis, qp, X);
350 else if (numGradDims_ == 3)
352 for (
int qp(0); qp < numQP; ++qp)
354 tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
355 vector_(cell, qp, X);
356 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
357 vector_(cell, qp, Y);
358 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
359 vector_(cell, qp, Z);
360 for (
int basis(0); basis < numBases; ++basis)
362 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
363 tmp[Z] * basis_(cell, basis, qp, Y);
364 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
365 tmp[X] * basis_(cell, basis, qp, Z);
366 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
367 tmp[Y] * basis_(cell, basis, qp, X);
374 const int numFieldMults(kokkosFieldMults_.extent(0));
375 if (numGradDims_ == 1)
377 for (
int qp(0); qp < numQP; ++qp)
379 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
380 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
381 for (
int fm(0); fm < numFieldMults; ++fm)
383 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
384 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
386 for (
int basis(0); basis < numBases; ++basis)
388 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
389 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
393 else if (numGradDims_ == 2)
395 for (
int qp(0); qp < numQP; ++qp)
397 tmp[X] = multiplier_ * vector_(cell, qp, X);
398 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
399 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
400 for (
int fm(0); fm < numFieldMults; ++fm)
402 tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
403 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
404 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
406 for (
int basis(0); basis < numBases; ++basis)
408 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
409 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
410 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
411 tmp[Y] * basis_(cell, basis, qp, X);
415 else if (numGradDims_ == 3)
417 for (
int qp(0); qp < numQP; ++qp)
419 tmp[X] = multiplier_ * vector_(cell, qp, X);
420 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
421 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
422 for (
int fm(0); fm < numFieldMults; ++fm)
424 tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
425 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
426 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
428 for (
int basis(0); basis < numBases; ++basis)
430 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
431 tmp[Z] * basis_(cell, basis, qp, Y);
432 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
433 tmp[X] * basis_(cell, basis, qp, Z);
434 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
435 tmp[Y] * basis_(cell, basis, qp, X);
447 template<
typename EvalT,
typename Traits>
453 using Kokkos::parallel_for;
454 using Kokkos::RangePolicy;
457 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
462 if (fieldMults_.size() == 0)
464 else if (fieldMults_.size() == 1)
475 template<
typename EvalT,
typename TRAITS>
476 Teuchos::RCP<Teuchos::ParameterList>
482 using PHX::DataLayout;
485 using Teuchos::ParameterList;
490 RCP<ParameterList> p = rcp(
new ParameterList);
492 RCP<const vector<string>> resNames;
493 p->set(
"Residual Names", resNames);
494 p->set<
string>(
"Vector Name",
"?");
495 RCP<BasisIRLayout> basis;
496 p->set(
"Basis", basis);
497 RCP<IntegrationRule> ir;
499 p->set<
double>(
"Multiplier", 1.0);
500 RCP<const vector<string>> fms;
501 p->set(
"Field Multipliers", fms);
502 RCP<DataLayout> vecDL;
503 p->set(
"Data Layout Vector", vecDL);
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.
typename EvalT::ScalarT ScalarT
The scalar type.
Kokkos::View< PHX::MDField< ScalarT, Cell, BASIS > * > fields_
The fields to which we'll contribute, or in which we'll store, the result of computing this integral.
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.
Kokkos::View< Kokkos::View< const ScalarT **, typename PHX::DevLayout< ScalarT >::type, PHX::Device > * > kokkosFieldMults_
The Kokkos::View representation of the (possibly empty) list of fields that are multipliers out in fr...
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.
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.
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.
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack... dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
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_