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,
78 const double& multiplier,
79 const std::vector<std::string>& fmNames,
81 const Teuchos::RCP<PHX::DataLayout>& vecDL )
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_);
141 fields_ =
OuterView(
"Integrator_GradBasisCrossVector::fields_", resNames.size());
144 for (
const auto& name : resNames)
148 for (std::size_t i=0; i<
fields_.extent(0); ++i) {
151 this->addContributedField(field);
153 this->addEvaluatedField(field);
159 kokkosFieldMults_ = PHX::View<PHX::UnmanagedView<const ScalarT**>*>(
"GradBasisCrossVector::KokkosFieldMultipliers", fmNames.size());
160 for (
const auto& name : fmNames)
167 string n(
"Integrator_GradBasisCrossVector (");
174 n += resNames[j] +
", ";
175 n += resNames[resNames.size()-1] +
"}";
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>
224 using Kokkos::createDynRankView;
228 auto fields_host_mirror_ = Kokkos::create_mirror_view(
fields_);
230 fields_host_mirror_(i) =
fields_host_[i].get_static_view();
232 Kokkos::deep_copy(
fields_,fields_host_mirror_);
237 field_mults_host_mirror_(i) =
fieldMults_[i].get_static_view();
249 template<
typename EvalT,
typename Traits>
250 template<
int NUM_FIELD_MULT>
251 KOKKOS_INLINE_FUNCTION
256 const size_t& cell)
const
261 const int numBases(
fields_[0].extent(1)), numQP(
vector_.extent(1));
263 for (
int dim(0); dim <
numDims_; ++dim)
264 for (
int basis(0); basis < numBases; ++basis)
265 fields_[dim](cell, basis) = 0.0;
270 const int X(0), Y(1), Z(2);
271 if (NUM_FIELD_MULT == 0)
275 for (
int qp(0); qp < numQP; ++qp)
279 for (
int basis(0); basis < numBases; ++basis)
281 fields_[Y](cell, basis) += tmp[Z] *
basis_(cell, basis, qp, X);
282 fields_[Z](cell, basis) += -tmp[Y] *
basis_(cell, basis, qp, X);
288 for (
int qp(0); qp < numQP; ++qp)
293 for (
int basis(0); basis < numBases; ++basis)
295 fields_[X](cell, basis) += -tmp[Z] *
basis_(cell, basis, qp, Y);
296 fields_[Y](cell, basis) += tmp[Z] *
basis_(cell, basis, qp, X);
297 fields_[Z](cell, basis) += tmp[X] *
basis_(cell, basis, qp, Y) -
298 tmp[Y] *
basis_(cell, basis, qp, X);
304 for (
int qp(0); qp < numQP; ++qp)
309 for (
int basis(0); basis < numBases; ++basis)
311 fields_[X](cell, basis) += tmp[Y] *
basis_(cell, basis, qp, Z) -
312 tmp[Z] *
basis_(cell, basis, qp, Y);
313 fields_[Y](cell, basis) += tmp[Z] *
basis_(cell, basis, qp, X) -
314 tmp[X] *
basis_(cell, basis, qp, Z);
315 fields_[Z](cell, basis) += tmp[X] *
basis_(cell, basis, qp, Y) -
316 tmp[Y] *
basis_(cell, basis, qp, X);
321 else if (NUM_FIELD_MULT == 1)
325 for (
int qp(0); qp < numQP; ++qp)
331 for (
int basis(0); basis < numBases; ++basis)
333 fields_[Y](cell, basis) += tmp[Z] *
basis_(cell, basis, qp, X);
334 fields_[Z](cell, basis) += -tmp[Y] *
basis_(cell, basis, qp, X);
340 for (
int qp(0); qp < numQP; ++qp)
348 for (
int basis(0); basis < numBases; ++basis)
350 fields_[X](cell, basis) += -tmp[Z] *
basis_(cell, basis, qp, Y);
351 fields_[Y](cell, basis) += tmp[Z] *
basis_(cell, basis, qp, X);
352 fields_[Z](cell, basis) += tmp[X] *
basis_(cell, basis, qp, Y) -
353 tmp[Y] *
basis_(cell, basis, qp, X);
359 for (
int qp(0); qp < numQP; ++qp)
367 for (
int basis(0); basis < numBases; ++basis)
369 fields_[X](cell, basis) += tmp[Y] *
basis_(cell, basis, qp, Z) -
370 tmp[Z] *
basis_(cell, basis, qp, Y);
371 fields_[Y](cell, basis) += tmp[Z] *
basis_(cell, basis, qp, X) -
372 tmp[X] *
basis_(cell, basis, qp, Z);
373 fields_[Z](cell, basis) += tmp[X] *
basis_(cell, basis, qp, Y) -
374 tmp[Y] *
basis_(cell, basis, qp, X);
384 for (
int qp(0); qp < numQP; ++qp)
388 for (
int fm(0); fm < numFieldMults; ++fm)
393 for (
int basis(0); basis < numBases; ++basis)
395 fields_[Y](cell, basis) += tmp[Z] *
basis_(cell, basis, qp, X);
396 fields_[Z](cell, basis) += -tmp[Y] *
basis_(cell, basis, qp, X);
402 for (
int qp(0); qp < numQP; ++qp)
407 for (
int fm(0); fm < numFieldMults; ++fm)
413 for (
int basis(0); basis < numBases; ++basis)
415 fields_[X](cell, basis) += -tmp[Z] *
basis_(cell, basis, qp, Y);
416 fields_[Y](cell, basis) += tmp[Z] *
basis_(cell, basis, qp, X);
417 fields_[Z](cell, basis) += tmp[X] *
basis_(cell, basis, qp, Y) -
418 tmp[Y] *
basis_(cell, basis, qp, X);
424 for (
int qp(0); qp < numQP; ++qp)
429 for (
int fm(0); fm < numFieldMults; ++fm)
435 for (
int basis(0); basis < numBases; ++basis)
437 fields_[X](cell, basis) += tmp[Y] *
basis_(cell, basis, qp, Z) -
438 tmp[Z] *
basis_(cell, basis, qp, Y);
439 fields_[Y](cell, basis) += tmp[Z] *
basis_(cell, basis, qp, X) -
440 tmp[X] *
basis_(cell, basis, qp, Z);
441 fields_[Z](cell, basis) += tmp[X] *
basis_(cell, basis, qp, Y) -
442 tmp[Y] *
basis_(cell, basis, qp, X);
454 template<
typename EvalT,
typename Traits>
460 using Kokkos::parallel_for;
461 using Kokkos::RangePolicy;
482 template<
typename EvalT,
typename TRAITS>
483 Teuchos::RCP<Teuchos::ParameterList>
489 using PHX::DataLayout;
492 using Teuchos::ParameterList;
497 RCP<ParameterList> p = rcp(
new ParameterList);
499 RCP<const vector<string>> resNames;
500 p->set(
"Residual Names", resNames);
501 p->set<
string>(
"Vector Name",
"?");
502 RCP<BasisIRLayout> basis;
503 p->set(
"Basis", basis);
504 RCP<IntegrationRule> ir;
506 p->set<
double>(
"Multiplier", 1.0);
507 RCP<const vector<string>> fms;
508 p->set(
"Field Multipliers", fms);
509 RCP<DataLayout> vecDL;
510 p->set(
"Data Layout Vector", vecDL);
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
WorksetDetailsAccessor wda
PHX::MDField< const ScalarT, Cell, IP, Dim > vector_
A field representing the vector-valued function we're integrating ( ).
PHX::View< InnerView * > OuterView
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...
std::size_t basisIndex_
The index in the Workset bases for our particular BasisIRLayout name.
typename EvalT::ScalarT ScalarT
The scalar type.
std::string basisName_
The name of the basis we're using.
PHX::MDField< double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dim > basis_
The gradient vector basis information necessary for integration.
double multiplier_
The scalar multiplier 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.
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_
const panzer::Workset & EvalData