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,
83 evalStyle_(evalStyle),
84 multiplier_(multiplier),
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;
106 "Integrator_GradBasisCrossVector called with the number of residual " \
107 "names not equal to three.")
108 for (
const auto& name : resNames)
110 "Integrator_GradBasisCrossVector called with an empty residual name.")
112 "Integrator_GradBasisCrossVector called with an empty vector name.")
115 "Error: Integrator_GradBasisCrossVector: Basis of type \""
116 << tmpBasis->name() <<
"\" does not support the gradient operator.")
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.")
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_);
141 for (
const auto& name : resNames)
143 MDField<ScalarT, Cell, BASIS> res(name, basis.
functional);
148 this->addContributedField(
field);
150 this->addEvaluatedField(
field);
156 typename DevLayout<ScalarT>::type, Device>*>(
157 "GradBasisCrossVector::KokkosFieldMultipliers", fmNames.size());
158 for (
const auto& name : fmNames)
165 string n(
"Integrator_GradBasisCrossVector (");
171 for (
size_t j=0; j < fields_.size() - 1; ++j)
172 n += fields_[j].fieldTag().name() +
", ";
173 n += fields_.back().fieldTag().name() +
"}";
182 template<
typename EvalT,
typename Traits>
189 p.get<
const std::vector<std::string>>(
"Residual Names"),
190 p.get<std::string>(
"Vector Name"),
193 p.get<double>(
"Multiplier"),
194 p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
195 (
"Field Multipliers") ?
196 (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
197 (
"Field Multipliers")) : std::vector<std::string>(),
198 p.isType<Teuchos::
RCP<PHX::DataLayout>>(
"Data Layout Vector") ?
199 p.get<Teuchos::
RCP<PHX::DataLayout>>(
"Data Layout Vector") :
215 template<
typename EvalT,
typename Traits>
227 for (
size_t i(0); i < fieldMults_.size(); ++i)
228 kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
240 template<
typename EvalT,
typename Traits>
241 template<
int NUM_FIELD_MULT>
242 KOKKOS_INLINE_FUNCTION
247 const size_t& cell)
const
252 const int numBases(fields_[0].extent(1)), numQP(vector_.extent(1));
254 for (
int dim(0); dim < numDims_; ++dim)
255 for (
int basis(0); basis < numBases; ++basis)
256 fields_[dim](cell, basis) = 0.0;
261 const int X(0), Y(1), Z(2);
262 if (NUM_FIELD_MULT == 0)
264 if (numGradDims_ == 1)
266 for (
int qp(0); qp < numQP; ++qp)
268 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
269 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
270 for (
int basis(0); basis < numBases; ++basis)
272 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
273 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
277 else if (numGradDims_ == 2)
279 for (
int qp(0); qp < numQP; ++qp)
281 tmp[X] = multiplier_ * vector_(cell, qp, X);
282 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
283 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
284 for (
int basis(0); basis < numBases; ++basis)
286 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
287 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
288 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
289 tmp[Y] * basis_(cell, basis, qp, X);
293 else if (numGradDims_ == 3)
295 for (
int qp(0); qp < numQP; ++qp)
297 tmp[X] = multiplier_ * vector_(cell, qp, X);
298 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
299 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
300 for (
int basis(0); basis < numBases; ++basis)
302 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
303 tmp[Z] * basis_(cell, basis, qp, Y);
304 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
305 tmp[X] * basis_(cell, basis, qp, Z);
306 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
307 tmp[Y] * basis_(cell, basis, qp, X);
312 else if (NUM_FIELD_MULT == 1)
314 if (numGradDims_ == 1)
316 for (
int qp(0); qp < numQP; ++qp)
318 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
319 vector_(cell, qp, Y);
320 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
321 vector_(cell, qp, Z);
322 for (
int basis(0); basis < numBases; ++basis)
324 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
325 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
329 else if (numGradDims_ == 2)
331 for (
int qp(0); qp < numQP; ++qp)
333 tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
334 vector_(cell, qp, X);
335 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
336 vector_(cell, qp, Y);
337 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
338 vector_(cell, qp, Z);
339 for (
int basis(0); basis < numBases; ++basis)
341 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
342 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
343 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
344 tmp[Y] * basis_(cell, basis, qp, X);
348 else if (numGradDims_ == 3)
350 for (
int qp(0); qp < numQP; ++qp)
352 tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
353 vector_(cell, qp, X);
354 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
355 vector_(cell, qp, Y);
356 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
357 vector_(cell, qp, Z);
358 for (
int basis(0); basis < numBases; ++basis)
360 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
361 tmp[Z] * basis_(cell, basis, qp, Y);
362 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
363 tmp[X] * basis_(cell, basis, qp, Z);
364 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
365 tmp[Y] * basis_(cell, basis, qp, X);
372 const int numFieldMults(kokkosFieldMults_.extent(0));
373 if (numGradDims_ == 1)
375 for (
int qp(0); qp < numQP; ++qp)
377 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
378 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
379 for (
int fm(0); fm < numFieldMults; ++fm)
381 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
382 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
384 for (
int basis(0); basis < numBases; ++basis)
386 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
387 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
391 else if (numGradDims_ == 2)
393 for (
int qp(0); qp < numQP; ++qp)
395 tmp[X] = multiplier_ * vector_(cell, qp, X);
396 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
397 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
398 for (
int fm(0); fm < numFieldMults; ++fm)
400 tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
401 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
402 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
404 for (
int basis(0); basis < numBases; ++basis)
406 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
407 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
408 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
409 tmp[Y] * basis_(cell, basis, qp, X);
413 else if (numGradDims_ == 3)
415 for (
int qp(0); qp < numQP; ++qp)
417 tmp[X] = multiplier_ * vector_(cell, qp, X);
418 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
419 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
420 for (
int fm(0); fm < numFieldMults; ++fm)
422 tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
423 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
424 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
426 for (
int basis(0); basis < numBases; ++basis)
428 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
429 tmp[Z] * basis_(cell, basis, qp, Y);
430 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
431 tmp[X] * basis_(cell, basis, qp, Z);
432 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
433 tmp[Y] * basis_(cell, basis, qp, X);
445 template<
typename EvalT,
typename Traits>
451 using Kokkos::parallel_for;
452 using Kokkos::RangePolicy;
455 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
460 if (fieldMults_.size() == 0)
462 else if (fieldMults_.size() == 1)
473 template<
typename EvalT,
typename TRAITS>
480 using PHX::DataLayout;
491 p->set(
"Residual Names", resNames);
492 p->set<
string>(
"Vector Name",
"?");
494 p->set(
"Basis", basis);
497 p->set<
double>(
"Multiplier", 1.0);
499 p->set(
"Field Multipliers", fms);
501 p->set(
"Data Layout Vector", vecDL);
508 #endif // PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP
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.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the 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.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
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...
EvaluatorStyle
An indication of how an Evaluator will behave.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
Teuchos::RCP< const PureBasis > getBasis() const
PHX::MDField< const ScalarT, Cell, IP, Dim > vector_
A field representing the vector-valued function we're integrating ( ).
double multiplier
The scalar multiplier out in front of the integral ( ).
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
typename EvalT::ScalarT ScalarT
The scalar type.
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.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
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...
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
std::vector< 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...
int numDims_
The number of dimensions associated with the vector.
int numGradDims_
The number of dimensions associated with the gradient.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_