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
76 evalStyle_(evalStyle),
77 multiplier_(multiplier),
78 basisName_(basis.name())
87 using std::invalid_argument;
88 using std::logic_error;
94 "Integrator_BasisTimesScalar called with an empty residual name.")
96 "Integrator_BasisTimesScalar called with an empty value name.")
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_);
118 View<Kokkos::View<const ScalarT**, typename PHX::DevLayout<ScalarT>::type, Kokkos::MemoryUnmanaged>*>(
"BasisTimesScalar::KokkosFieldMultipliers",
120 for (
const auto& name : fmNames)
127 string n(
"Integrator_BasisTimesScalar (");
132 n +=
"): " +
field_.fieldTag().name();
141 template<
typename EvalT,
typename Traits>
148 p.get<std::string>(
"Residual Name"),
149 p.get<std::string>(
"Value Name"),
152 p.get<double>(
"Multiplier"),
153 p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
154 (
"Field Multipliers") ?
155 (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
156 (
"Field Multipliers")) : std::vector<std::string>())
171 template<
typename EvalT,
typename Traits>
181 auto kokkosFieldMults_h = Kokkos::create_mirror_view(kokkosFieldMults_);
184 for (
size_t i(0); i < fieldMults_.size(); ++i)
185 kokkosFieldMults_h(i) = fieldMults_[i].get_static_view();
187 Kokkos::deep_copy(kokkosFieldMults_, kokkosFieldMults_h);
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);
196 if (fieldMults_.size() > 1)
197 tmp2_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",field_.extent(0),fadSize);
199 tmp_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",field_.extent(0));
200 if (fieldMults_.size() > 1)
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)
236 tmp_(cell) = multiplier_ * scalar_(cell, qp);
237 for (
int basis(0); basis < numBases; ++basis)
238 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
241 else if (NUM_FIELD_MULT == 1)
246 for (
int qp(0); qp < numQP; ++qp)
248 tmp_(cell) = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
249 for (
int basis(0); basis < numBases; ++basis)
250 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
259 const int numFieldMults(kokkosFieldMults_.extent(0));
260 for (
int qp(0); qp < numQP; ++qp)
263 for (
int fm(0); fm < numFieldMults; ++fm)
264 tmp2_(cell) *= kokkosFieldMults_(fm)(cell, qp);
265 tmp_(cell) = multiplier_ * scalar_(cell, qp) * tmp2_(cell);
266 for (
int basis(0); basis < numBases; ++basis)
267 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
277 template<
typename EvalT,
typename Traits>
283 using Kokkos::parallel_for;
284 using Kokkos::RangePolicy;
287 basis_ = this->wda(workset).bases[basisIndex_]->weighted_basis_scalar;
292 if (fieldMults_.size() == 0)
294 else if (fieldMults_.size() == 1)
305 template<
typename EvalT,
typename TRAITS>
320 p->set<
string>(
"Residual Name",
"?");
321 p->set<
string>(
"Value Name",
"?");
323 p->set(
"Basis", basis);
326 p->set<
double>(
"Multiplier", 1.0);
328 p->set(
"Field Multipliers", fms);
334 #endif // PANZER_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP
int num_cells
DEPRECATED - use: numCells()
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...
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.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
Teuchos::RCP< const PureBasis > getBasis() 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.
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.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > scalar_
A field representing the scalar function we're integrating ( ).
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
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...
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
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< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_