11 #ifndef PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
12 #define PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
21 #include "Intrepid2_FunctionSpaceTools.hpp"
24 #include "Kokkos_ViewFactory.hpp"
38 template<
typename EvalT,
typename Traits>
42 const std::vector<std::string>& resNames,
43 const std::string& scalar,
47 const std::vector<std::string>& fmNames)
50 evalStyle_(evalStyle),
51 multiplier_(multiplier),
52 numDims_(ir.dl_vector->extent(2)),
53 basisName_(basis.name())
64 using std::invalid_argument;
65 using std::logic_error;
70 for (
const auto& name : resNames)
72 "Integrator_GradBasisTimesScalar called with an empty residual name.")
74 "Integrator_GradBasisTimesScalar called with an empty scalar name.")
76 logic_error,
"Error: Integrator_GradBasisTimesScalar: The number of " \
77 "residuals must equal the number of gradient components (mesh " \
81 "Error: Integrator_GradBasisTimesScalar: Basis of type \""
82 << tmpBasis->name() <<
"\" does not support the gradient operation.")
86 this->addDependentField(
scalar_);
91 fields_ =
OuterView(
"Integrator_GradBasisCrossVector::fields_", resNames.size());
94 for (
const auto& name : resNames)
97 for (std::size_t i=0; i<
fields_.extent(0); ++i) {
100 this->addContributedField(
field);
102 this->addEvaluatedField(
field);
109 "GradBasisTimesScalar::KokkosFieldMultipliers", fmNames.size());
110 for (
const auto& name : fmNames)
117 string n(
"Integrator_GradBasisTimesScalar (");
124 n += resNames[j] +
", ";
125 n += resNames[resNames.size()-1] +
"}";
134 template<
typename EvalT,
typename Traits>
141 p.get<
const std::vector<std::string>>(
"Residual Names"),
142 p.get<std::string>(
"Scalar Name"),
145 p.get<double>(
"Multiplier"),
146 p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
147 (
"Field Multipliers") ?
148 (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
149 (
"Field Multipliers")) : std::vector<std::string>())
164 template<
typename EvalT,
typename Traits>
175 auto fields_host_mirror = Kokkos::create_mirror_view(fields_);
176 for (
size_t i(0); i < fields_host_.size(); ++i)
177 fields_host_mirror(i) = fields_host_[i].get_static_view();
178 Kokkos::deep_copy(fields_,fields_host_mirror);
181 auto field_mults_host_mirror = Kokkos::create_mirror_view(kokkosFieldMults_);
182 for (
size_t i(0); i < fieldMults_.size(); ++i)
183 field_mults_host_mirror(i) = fieldMults_[i].get_static_view();
184 Kokkos::deep_copy(kokkosFieldMults_,field_mults_host_mirror);
190 template<
typename EvalT,
typename Traits>
191 template<
int NUM_FIELD_MULT>
192 KOKKOS_INLINE_FUNCTION
197 const size_t& cell)
const
202 const int numQP(scalar_.extent(1)), numBases(fields_[0].extent(1));
204 for (
int dim(0); dim < numDims_; ++dim)
205 for (
int basis(0); basis < numBases; ++basis)
206 fields_[dim](cell, basis) = 0.0;
211 if (NUM_FIELD_MULT == 0)
216 for (
int qp(0); qp < numQP; ++qp)
218 tmp = multiplier_ * scalar_(cell, qp);
219 for (
int basis(0); basis < numBases; ++basis)
220 for (
int dim(0); dim < numDims_; ++dim)
221 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
224 else if (NUM_FIELD_MULT == 1)
229 for (
int qp(0); qp < numQP; ++qp)
231 tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
232 for (
int basis(0); basis < numBases; ++basis)
233 for (
int dim(0); dim < numDims_; ++dim)
234 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
243 const int numFieldMults(kokkosFieldMults_.extent(0));
244 for (
int qp(0); qp < numQP; ++qp)
247 for (
int fm(0); fm < numFieldMults; ++fm)
248 fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
249 tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
250 for (
int basis(0); basis < numBases; ++basis)
251 for (
int dim(0); dim < numDims_; ++dim)
252 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
262 template<
typename EvalT,
typename Traits>
268 using Kokkos::parallel_for;
269 using Kokkos::RangePolicy;
272 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
277 if (fieldMults_.size() == 0)
279 else if (fieldMults_.size() == 1)
290 template<
typename EvalT,
typename TRAITS>
307 p->set(
"Residual Names", resNames);
308 p->set<
string>(
"Scalar Name",
"?");
310 p->set(
"Basis", basis);
313 p->set<
double>(
"Multiplier", 1.0);
315 p->set(
"Field Multipliers", fms);
322 #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.
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 ( ).
int num_cells
DEPRECATED - use: numCells()
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
EvaluatorStyle
An indication of how an Evaluator will behave.
PHX::View< InnerView * > OuterView
Teuchos::RCP< const PureBasis > getBasis() const
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
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.
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 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_