43 #ifndef PANZER_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP 44 #define PANZER_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP 64 template<
typename EvalT,
typename Traits>
68 const std::string& resName,
69 const std::string& valName,
73 const std::vector<std::string>& fmNames
78 basisName_(basis.name())
87 using PHX::typeAsString;
88 using std::invalid_argument;
89 using std::logic_error;
95 "Integrator_BasisTimesScalar called with an empty residual name.")
97 "Integrator_BasisTimesScalar called with an empty value name.")
100 "Error: Integrator_BasisTimesScalar: Basis of type \"" 101 << tmpBasis->name() <<
"\" is not a scalar basis.")
105 this->addDependentField(
scalar_);
111 this->addContributedField(
field_);
113 this->addEvaluatedField(
field_);
119 View<View<const ScalarT**>*>(
"BasisTimesScalar::KokkosFieldMultipliers",
121 for (
const auto& name : fmNames)
128 string n(
"Integrator_BasisTimesScalar (");
133 n +=
"): " +
field_.fieldTag().name();
142 template<
typename EvalT,
typename Traits>
149 p.get<
std::string>(
"Residual Name"),
150 p.get<
std::string>(
"Value Name"),
153 p.get<double>(
"Multiplier"),
155 (
"Field Multipliers") ?
157 (
"Field Multipliers")) :
std::vector<
std::string>())
172 template<
typename EvalT,
typename Traits>
183 for (
size_t i(0); i < fieldMults_.size(); ++i)
184 kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
187 numQP_ = scalar_.extent(1);
198 template<
typename EvalT,
typename Traits>
199 template<
int NUM_FIELD_MULT>
200 KOKKOS_INLINE_FUNCTION
205 const size_t& cell)
const 210 const int numBases(basis_.extent(1));
212 for (
int basis(0); basis < numBases; ++basis)
213 field_(cell, basis) = 0.0;
218 if (NUM_FIELD_MULT == 0)
223 for (
int qp(0); qp < numQP_; ++qp)
225 tmp = multiplier_ * scalar_(cell, qp);
226 for (
int basis(0); basis < numBases; ++basis)
227 field_(cell, basis) += basis_(cell, basis, qp) * tmp;
230 else if (NUM_FIELD_MULT == 1)
235 for (
int qp(0); qp < numQP_; ++qp)
237 tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
238 for (
int basis(0); basis < numBases; ++basis)
239 field_(cell, basis) += basis_(cell, basis, qp) * tmp;
248 const int numFieldMults(kokkosFieldMults_.extent(0));
249 for (
int qp(0); qp < numQP_; ++qp)
252 for (
int fm(0); fm < numFieldMults; ++fm)
253 fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
254 tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
255 for (
int basis(0); basis < numBases; ++basis)
256 field_(cell, basis) += basis_(cell, basis, qp) * tmp;
266 template<
typename EvalT,
typename Traits>
272 using Kokkos::parallel_for;
273 using Kokkos::RangePolicy;
276 basis_ = this->wda(workset).bases[basisIndex_]->weighted_basis_scalar;
281 if (fieldMults_.size() == 0)
283 else if (fieldMults_.size() == 1)
294 template<
typename EvalT,
typename TRAITS>
309 p->set<
string>(
"Residual Name",
"?");
310 p->set<
string>(
"Value Name",
"?");
312 p->set(
"Basis", basis);
315 p->set<
double>(
"Multiplier", 1.0);
317 p->set(
"Field Multipliers", fms);
323 #endif // PANZER_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
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...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Teuchos::RCP< const PureBasis > getBasis() const
Integrator_BasisTimesScalar(const panzer::EvaluatorStyle &evalStyle, const std::string &resName, const std::string &valName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >())
Main Constructor.
EvaluatorStyle
An indication of how an Evaluator will behave.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
double multiplier
The scalar multiplier out in front of the integral ( ).
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
Kokkos::View< Kokkos::View< const ScalarT ** > * > kokkosFieldMults_
The Kokkos::View representation of the (possibly empty) list of fields that are multipliers out in fr...
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
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.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > scalar_
A field representing the scalar function we're integrating ( ).
typename EvalT::ScalarT ScalarT
The scalar type.
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Description and data layouts associated with a particular basis.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ...
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_