11#ifndef PANZER_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP
12#define PANZER_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP
32 template<
typename EvalT,
typename Traits>
36 const std::string& resName,
37 const std::string& valName,
41 const std::vector<std::string>& fmNames
46 basisName_(basis.name())
55 using std::invalid_argument;
56 using std::logic_error;
61 TEUCHOS_TEST_FOR_EXCEPTION(resName ==
"", invalid_argument,
"Error: " \
62 "Integrator_BasisTimesScalar called with an empty residual name.")
63 TEUCHOS_TEST_FOR_EXCEPTION(valName ==
"", invalid_argument,
"Error: " \
64 "Integrator_BasisTimesScalar called with an empty value name.")
65 RCP<const PureBasis> tmpBasis = basis.
getBasis();
66 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isScalarBasis(), logic_error,
67 "Error: Integrator_BasisTimesScalar: Basis of type \""
68 << tmpBasis->name() <<
"\" is not a scalar basis.")
72 this->addDependentField(
scalar_);
78 this->addContributedField(
field_);
80 this->addEvaluatedField(
field_);
86 View<Kokkos::View<const ScalarT**, typename PHX::DevLayout<ScalarT>::type, Kokkos::MemoryUnmanaged>*>(
"BasisTimesScalar::KokkosFieldMultipliers",
88 for (
const auto& name : fmNames)
95 string n(
"Integrator_BasisTimesScalar (");
100 n +=
"): " +
field_.fieldTag().name();
109 template<
typename EvalT,
typename Traits>
115 p.get<
std::string>(
"Residual Name"),
116 p.get<
std::string>(
"Value Name"),
119 p.get<double>(
"Multiplier"),
121 (
"Field Multipliers") ?
123 (
"Field Multipliers")) :
std::vector<
std::string>())
125 using Teuchos::ParameterList;
130 p.validateParameters(*validParams);
138 template<
typename EvalT,
typename Traits>
148 auto kokkosFieldMults_h = Kokkos::create_mirror_view(kokkosFieldMults_);
151 for (
size_t i(0); i < fieldMults_.size(); ++i)
152 kokkosFieldMults_h(i) = fieldMults_[i].get_static_view();
154 Kokkos::deep_copy(kokkosFieldMults_, kokkosFieldMults_h);
160 if (Sacado::IsADType<ScalarT>::value) {
161 const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
162 tmp_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",field_.extent(0),fadSize);
163 if (fieldMults_.size() > 1)
164 tmp2_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",field_.extent(0),fadSize);
166 tmp_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",field_.extent(0));
167 if (fieldMults_.size() > 1)
168 tmp2_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",field_.extent(0));
177 template<
typename EvalT,
typename Traits>
178 template<
int NUM_FIELD_MULT>
179 KOKKOS_INLINE_FUNCTION
184 const size_t& cell)
const
189 const int numQP(scalar_.extent(1)), numBases(basis_.extent(1));
191 for (
int basis(0); basis < numBases; ++basis)
192 field_(cell, basis) = 0.0;
196 if (NUM_FIELD_MULT == 0)
201 for (
int qp(0); qp < numQP; ++qp)
203 tmp_(cell) = multiplier_ * scalar_(cell, qp);
204 for (
int basis(0); basis < numBases; ++basis)
205 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
208 else if (NUM_FIELD_MULT == 1)
213 for (
int qp(0); qp < numQP; ++qp)
215 tmp_(cell) = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
216 for (
int basis(0); basis < numBases; ++basis)
217 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
226 const int numFieldMults(kokkosFieldMults_.extent(0));
227 for (
int qp(0); qp < numQP; ++qp)
230 for (
int fm(0); fm < numFieldMults; ++fm)
231 tmp2_(cell) *= kokkosFieldMults_(fm)(cell, qp);
232 tmp_(cell) = multiplier_ * scalar_(cell, qp) * tmp2_(cell);
233 for (
int basis(0); basis < numBases; ++basis)
234 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
244 template<
typename EvalT,
typename Traits>
250 using Kokkos::parallel_for;
251 using Kokkos::RangePolicy;
254 basis_ = this->wda(workset).bases[basisIndex_]->weighted_basis_scalar;
259 if (fieldMults_.size() == 0)
261 else if (fieldMults_.size() == 1)
272 template<
typename EvalT,
typename TRAITS>
273 Teuchos::RCP<Teuchos::ParameterList>
281 using Teuchos::ParameterList;
286 RCP<ParameterList> p = rcp(
new ParameterList);
287 p->set<
string>(
"Residual Name",
"?");
288 p->set<
string>(
"Value Name",
"?");
289 RCP<BasisIRLayout> basis;
290 p->set(
"Basis", basis);
291 RCP<IntegrationRule> ir;
293 p->set<
double>(
"Multiplier", 1.0);
294 RCP<const vector<string>> fms;
295 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 ( ).
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
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 ( ,...
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::MDField< const ScalarT, panzer::Cell, panzer::IP > scalar_
A field representing the scalar function we're integrating ( ).
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
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_