43 #ifndef PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP 44 #define PANZER_EVALUATOR_GRADBASISTIMESSCALAR_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& scalar,
79 const std::vector<std::string>& fmNames)
84 numDims_(ir.dl_vector->extent(2)),
85 basisName_(basis.name())
96 using std::invalid_argument;
97 using std::logic_error;
102 for (
const auto& name : resNames)
103 TEUCHOS_TEST_FOR_EXCEPTION(name ==
"", invalid_argument,
"Error: " \
104 "Integrator_GradBasisTimesScalar called with an empty residual name.")
105 TEUCHOS_TEST_FOR_EXCEPTION(scalar ==
"", invalid_argument,
"Error: " \
106 "Integrator_GradBasisTimesScalar called with an empty scalar name.")
107 TEUCHOS_TEST_FOR_EXCEPTION(
numDims_ != static_cast<int>(resNames.size()),
108 logic_error,
"Error: Integrator_GradBasisTimesScalar: The number of " \
109 "residuals must equal the number of gradient components (mesh " \
111 RCP<const PureBasis> tmpBasis = basis.
getBasis();
112 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsGrad(), logic_error,
113 "Error: Integrator_GradBasisTimesScalar: Basis of type \"" 114 << tmpBasis->name() <<
"\" does not support the gradient operation.")
118 this->addDependentField(
scalar_);
122 fields_ = PHX::View<PHX::MDField<ScalarT, Cell, BASIS>*>(
"Integrator_GradBasisTimesScalar",resNames.size());
125 for (
const auto& name : resNames)
128 for (std::size_t i=0; i<
fields_.extent(0); ++i) {
131 this->addContributedField(
field);
133 this->addEvaluatedField(
field);
140 "GradBasisTimesScalar::KokkosFieldMultipliers", fmNames.size());
141 for (
const auto& name : fmNames)
148 string n(
"Integrator_GradBasisTimesScalar (");
154 for (
size_t j=0; j <
fields_.extent(0) - 1; ++j)
155 n +=
fields_[j].fieldTag().name() +
", ";
165 template<
typename EvalT,
typename Traits>
168 const Teuchos::ParameterList& p)
172 p.get<const
std::vector<
std::string>>(
"Residual Names"),
173 p.get<
std::string>(
"Scalar Name"),
176 p.get<double>(
"Multiplier"),
178 (
"Field Multipliers") ?
180 (
"Field Multipliers")) :
std::vector<
std::string>())
182 using Teuchos::ParameterList;
187 p.validateParameters(*validParams);
195 template<
typename EvalT,
typename Traits>
206 for (
size_t i(0); i < fieldMults_.size(); ++i)
207 kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
213 template<
typename EvalT,
typename Traits>
214 template<
int NUM_FIELD_MULT>
215 KOKKOS_INLINE_FUNCTION
220 const size_t& cell)
const 225 const int numQP(scalar_.extent(1)), numBases(fields_[0].extent(1));
227 for (
int dim(0); dim < numDims_; ++dim)
228 for (
int basis(0); basis < numBases; ++basis)
229 fields_[dim](cell, basis) = 0.0;
234 if (NUM_FIELD_MULT == 0)
239 for (
int qp(0); qp < numQP; ++qp)
241 tmp = multiplier_ * scalar_(cell, qp);
242 for (
int basis(0); basis < numBases; ++basis)
243 for (
int dim(0); dim < numDims_; ++dim)
244 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
247 else if (NUM_FIELD_MULT == 1)
252 for (
int qp(0); qp < numQP; ++qp)
254 tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
255 for (
int basis(0); basis < numBases; ++basis)
256 for (
int dim(0); dim < numDims_; ++dim)
257 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
266 const int numFieldMults(kokkosFieldMults_.extent(0));
267 for (
int qp(0); qp < numQP; ++qp)
270 for (
int fm(0); fm < numFieldMults; ++fm)
271 fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
272 tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
273 for (
int basis(0); basis < numBases; ++basis)
274 for (
int dim(0); dim < numDims_; ++dim)
275 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
285 template<
typename EvalT,
typename Traits>
291 using Kokkos::parallel_for;
292 using Kokkos::RangePolicy;
295 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
300 if (fieldMults_.size() == 0)
302 else if (fieldMults_.size() == 1)
313 template<
typename EvalT,
typename TRAITS>
314 Teuchos::RCP<Teuchos::ParameterList>
322 using Teuchos::ParameterList;
327 RCP<ParameterList> p = rcp(
new ParameterList);
329 RCP<const vector<string>> resNames;
330 p->set(
"Residual Names", resNames);
331 p->set<
string>(
"Scalar Name",
"?");
332 RCP<BasisIRLayout> basis;
333 p->set(
"Basis", basis);
334 RCP<IntegrationRule> ir;
336 p->set<
double>(
"Multiplier", 1.0);
337 RCP<const vector<string>> fms;
338 p->set(
"Field Multipliers", fms);
345 #endif // PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP typename EvalT::ScalarT ScalarT
The scalar type.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we're integrating ( ).
PHX::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...
int num_cells
DEPRECATED - use: numCells()
PHX::View< PHX::View< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
Teuchos::RCP< const PureBasis > getBasis() const
EvaluatorStyle
An indication of how an Evaluator will behave.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
double multiplier
The scalar multiplier out in front of the integral ( ).
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.
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The scalar multiplier out in front of the integral ( ).
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
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.
int numDims_
The number of dimensions associated with the gradient.
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...
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
Description and data layouts associated with a particular basis.
Integrator_GradBasisTimesScalar(const panzer::EvaluatorStyle &evalStyle, const std::vector< std::string > &resNames, const std::string &scalar, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >())
Main Constructor.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_