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,
72 const double& multiplier,
73 const std::vector<std::string>& fmNames
87 using std::invalid_argument;
88 using std::logic_error;
93 TEUCHOS_TEST_FOR_EXCEPTION(resName ==
"", invalid_argument,
"Error: " \
94 "Integrator_BasisTimesScalar called with an empty residual name.")
95 TEUCHOS_TEST_FOR_EXCEPTION(valName ==
"", invalid_argument,
"Error: " \
96 "Integrator_BasisTimesScalar called with an empty value name.")
97 RCP<const PureBasis> tmpBasis = basis.
getBasis();
98 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isScalarBasis(), logic_error,
99 "Error: Integrator_BasisTimesScalar: Basis of type \""
100 << tmpBasis->name() <<
"\" is not a scalar basis.")
104 this->addDependentField(
scalar_);
110 this->addContributedField(
field_);
112 this->addEvaluatedField(
field_);
120 for (
const auto& name : fmNames)
127 string n(
"Integrator_BasisTimesScalar (");
132 n +=
"): " +
field_.fieldTag().name();
141 template<
typename EvalT,
typename Traits>
144 const Teuchos::ParameterList& p)
148 p.get<
std::string>(
"Residual Name"),
149 p.get<
std::string>(
"Value Name"),
152 p.get<double>(
"Multiplier"),
154 (
"Field Multipliers") ?
156 (
"Field Multipliers")) :
std::vector<
std::string>())
158 using Teuchos::ParameterList;
163 p.validateParameters(*validParams);
171 template<
typename EvalT,
typename Traits>
185 kokkosFieldMults_h(i) =
fieldMults_[i].get_static_view();
193 if (Sacado::IsADType<ScalarT>::value) {
194 const auto fadSize = Kokkos::dimension_scalar(
field_.get_view());
195 tmp_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",
field_.extent(0),fadSize);
197 tmp2_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",
field_.extent(0),fadSize);
199 tmp_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",
field_.extent(0));
201 tmp2_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",
field_.extent(0));
210 template<
typename EvalT,
typename Traits>
211 template<
int NUM_FIELD_MULT>
212 KOKKOS_INLINE_FUNCTION
217 const size_t& cell)
const
222 const int numQP(
scalar_.extent(1)), numBases(
basis_.extent(1));
224 for (
int basis(0); basis < numBases; ++basis)
225 field_(cell, basis) = 0.0;
229 if (NUM_FIELD_MULT == 0)
234 for (
int qp(0); qp < numQP; ++qp)
237 for (
int basis(0); basis < numBases; ++basis)
241 else if (NUM_FIELD_MULT == 1)
246 for (
int qp(0); qp < numQP; ++qp)
249 for (
int basis(0); basis < numBases; ++basis)
260 for (
int qp(0); qp < numQP; ++qp)
263 for (
int fm(0); fm < numFieldMults; ++fm)
266 for (
int basis(0); basis < numBases; ++basis)
277 template<
typename EvalT,
typename Traits>
283 using Kokkos::parallel_for;
284 using Kokkos::RangePolicy;
305 template<
typename EvalT,
typename TRAITS>
306 Teuchos::RCP<Teuchos::ParameterList>
314 using Teuchos::ParameterList;
319 RCP<ParameterList> p = rcp(
new ParameterList);
320 p->set<
string>(
"Residual Name",
"?");
321 p->set<
string>(
"Value Name",
"?");
322 RCP<BasisIRLayout> basis;
323 p->set(
"Basis", basis);
324 RCP<IntegrationRule> ir;
326 p->set<
double>(
"Multiplier", 1.0);
327 RCP<const vector<string>> fms;
328 p->set(
"Field Multipliers", fms);
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
WorksetDetailsAccessor wda
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
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 ( ,...
PHX::View< ScalarT * > tmp2_
For storing temporaries, one value per thread.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
PHX::View< Kokkos::View< const ScalarT **, typename PHX::DevLayout< ScalarT >::type, Kokkos::MemoryUnmanaged > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
PHX::View< ScalarT * > tmp_
For storing temporaries, one value per thread.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > scalar_
A field representing the scalar function we're integrating ( ).
PHX::MDField< double, panzer::Cell, panzer::BASIS, panzer::IP > basis_
The scalar basis information necessary for integration.
double multiplier_
The scalar multiplier out in front of the integral ( ).
std::size_t basisIndex_
The index in the Workset bases for our particular BasisIRLayout name.
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
std::string basisName_
The name of the basis we're using.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
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.
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 > 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_
const panzer::Workset & EvalData