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_);
123 fields_ =
OuterView(
"Integrator_GradBasisCrossVector::fields_", resNames.size());
126 for (
const auto& name : resNames)
129 for (std::size_t i=0; i<
fields_.extent(0); ++i) {
132 this->addContributedField(
field);
134 this->addEvaluatedField(
field);
141 "GradBasisTimesScalar::KokkosFieldMultipliers", fmNames.size());
142 for (
const auto& name : fmNames)
149 string n(
"Integrator_GradBasisTimesScalar (");
156 n += resNames[j] +
", ";
157 n += resNames[resNames.size()-1] +
"}";
166 template<
typename EvalT,
typename Traits>
169 const Teuchos::ParameterList& p)
173 p.get<const
std::vector<
std::string>>(
"Residual Names"),
174 p.get<
std::string>(
"Scalar Name"),
177 p.get<double>(
"Multiplier"),
179 (
"Field Multipliers") ?
181 (
"Field Multipliers")) :
std::vector<
std::string>())
183 using Teuchos::ParameterList;
188 p.validateParameters(*validParams);
196 template<
typename EvalT,
typename Traits>
203 using Kokkos::createDynRankView;
207 auto fields_host_mirror = Kokkos::create_mirror_view(fields_);
208 for (
size_t i(0); i < fields_host_.size(); ++i)
209 fields_host_mirror(i) = fields_host_[i].get_static_view();
210 Kokkos::deep_copy(fields_,fields_host_mirror);
213 auto field_mults_host_mirror = Kokkos::create_mirror_view(kokkosFieldMults_);
214 for (
size_t i(0); i < fieldMults_.size(); ++i)
215 field_mults_host_mirror(i) = fieldMults_[i].get_static_view();
216 Kokkos::deep_copy(kokkosFieldMults_,field_mults_host_mirror);
222 template<
typename EvalT,
typename Traits>
223 template<
int NUM_FIELD_MULT>
224 KOKKOS_INLINE_FUNCTION
229 const size_t& cell)
const
234 const int numQP(scalar_.extent(1)), numBases(fields_[0].extent(1));
236 for (
int dim(0); dim < numDims_; ++dim)
237 for (
int basis(0); basis < numBases; ++basis)
238 fields_[dim](cell, basis) = 0.0;
243 if (NUM_FIELD_MULT == 0)
248 for (
int qp(0); qp < numQP; ++qp)
250 tmp = multiplier_ * scalar_(cell, qp);
251 for (
int basis(0); basis < numBases; ++basis)
252 for (
int dim(0); dim < numDims_; ++dim)
253 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
256 else if (NUM_FIELD_MULT == 1)
261 for (
int qp(0); qp < numQP; ++qp)
263 tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
264 for (
int basis(0); basis < numBases; ++basis)
265 for (
int dim(0); dim < numDims_; ++dim)
266 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
275 const int numFieldMults(kokkosFieldMults_.extent(0));
276 for (
int qp(0); qp < numQP; ++qp)
279 for (
int fm(0); fm < numFieldMults; ++fm)
280 fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
281 tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
282 for (
int basis(0); basis < numBases; ++basis)
283 for (
int dim(0); dim < numDims_; ++dim)
284 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
294 template<
typename EvalT,
typename Traits>
300 using Kokkos::parallel_for;
301 using Kokkos::RangePolicy;
304 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
309 if (fieldMults_.size() == 0)
311 else if (fieldMults_.size() == 1)
322 template<
typename EvalT,
typename TRAITS>
323 Teuchos::RCP<Teuchos::ParameterList>
331 using Teuchos::ParameterList;
336 RCP<ParameterList> p = rcp(
new ParameterList);
338 RCP<const vector<string>> resNames;
339 p->set(
"Residual Names", resNames);
340 p->set<
string>(
"Scalar Name",
"?");
341 RCP<BasisIRLayout> basis;
342 p->set(
"Basis", basis);
343 RCP<IntegrationRule> ir;
345 p->set<
double>(
"Multiplier", 1.0);
346 RCP<const vector<string>> fms;
347 p->set(
"Field Multipliers", fms);
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
double multiplier
The scalar multiplier out in front of the integral ( ).
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.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
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.
PHX::View< InnerView * > OuterView
int numDims_
The number of dimensions associated with the gradient.
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we're integrating ( ).
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.
typename EvalT::ScalarT ScalarT
The scalar type.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The scalar multiplier out in front of the integral ( ).
PHX::View< PHX::UnmanagedView< 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.
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
Description and data layouts associated with a particular basis.
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_