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)
82 evalStyle_(evalStyle),
83 multiplier_(multiplier),
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)
104 "Integrator_GradBasisTimesScalar called with an empty residual name.")
106 "Integrator_GradBasisTimesScalar called with an empty scalar name.")
108 logic_error,
"Error: Integrator_GradBasisTimesScalar: The number of " \
109 "residuals must equal the number of gradient components (mesh " \
113 "Error: Integrator_GradBasisTimesScalar: Basis of type \""
114 << tmpBasis->name() <<
"\" does not support the gradient operation.")
118 this->addDependentField(
scalar_);
122 fields_ = Kokkos::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 typename DevLayout<ScalarT>::type, Device>*>(
141 "GradBasisTimesScalar::KokkosFieldMultipliers", fmNames.size());
142 for (
const auto& name : fmNames)
149 string n(
"Integrator_GradBasisTimesScalar (");
155 for (
size_t j=0; j <
fields_.extent(0) - 1; ++j)
156 n +=
fields_[j].fieldTag().name() +
", ";
166 template<
typename EvalT,
typename Traits>
173 p.get<
const std::vector<std::string>>(
"Residual Names"),
174 p.get<std::string>(
"Scalar Name"),
177 p.get<double>(
"Multiplier"),
178 p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
179 (
"Field Multipliers") ?
180 (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
181 (
"Field Multipliers")) : std::vector<std::string>())
196 template<
typename EvalT,
typename Traits>
208 for (
size_t i(0); i < fieldMults_.size(); ++i)
209 kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
216 template<
typename EvalT,
typename Traits>
217 template<
int NUM_FIELD_MULT>
218 KOKKOS_INLINE_FUNCTION
223 const size_t& cell)
const
228 const int numQP(scalar_.extent(1)), numBases(fields_[0].extent(1));
230 for (
int dim(0); dim < numDims_; ++dim)
231 for (
int basis(0); basis < numBases; ++basis)
232 fields_[dim](cell, basis) = 0.0;
237 if (NUM_FIELD_MULT == 0)
242 for (
int qp(0); qp < numQP; ++qp)
244 tmp = multiplier_ * scalar_(cell, qp);
245 for (
int basis(0); basis < numBases; ++basis)
246 for (
int dim(0); dim < numDims_; ++dim)
247 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
250 else if (NUM_FIELD_MULT == 1)
255 for (
int qp(0); qp < numQP; ++qp)
257 tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
258 for (
int basis(0); basis < numBases; ++basis)
259 for (
int dim(0); dim < numDims_; ++dim)
260 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
269 const int numFieldMults(kokkosFieldMults_.extent(0));
270 for (
int qp(0); qp < numQP; ++qp)
273 for (
int fm(0); fm < numFieldMults; ++fm)
274 fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
275 tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
276 for (
int basis(0); basis < numBases; ++basis)
277 for (
int dim(0); dim < numDims_; ++dim)
278 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
288 template<
typename EvalT,
typename Traits>
294 using Kokkos::parallel_for;
295 using Kokkos::RangePolicy;
298 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
303 if (fieldMults_.size() == 0)
305 else if (fieldMults_.size() == 1)
316 template<
typename EvalT,
typename TRAITS>
333 p->set(
"Residual Names", resNames);
334 p->set<
string>(
"Scalar Name",
"?");
336 p->set(
"Basis", basis);
339 p->set<
double>(
"Multiplier", 1.0);
341 p->set(
"Field Multipliers", fms);
348 #endif // PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
typename EvalT::ScalarT ScalarT
The scalar type.
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...
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
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.
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we're integrating ( ).
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
EvaluatorStyle
An indication of how an Evaluator will behave.
Teuchos::RCP< const PureBasis > getBasis() const
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
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...
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
double multiplier
The scalar multiplier out in front of the integral ( ).
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
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...
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
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>
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_