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
44 evalStyle_(evalStyle),
45 multiplier_(multiplier),
46 basisName_(basis.name())
55 using std::invalid_argument;
56 using std::logic_error;
62 "Integrator_BasisTimesScalar called with an empty residual name.")
64 "Integrator_BasisTimesScalar called with an empty value name.")
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>
116 p.get<std::string>(
"Residual Name"),
117 p.get<std::string>(
"Value Name"),
120 p.get<double>(
"Multiplier"),
121 p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
122 (
"Field Multipliers") ?
123 (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
124 (
"Field Multipliers")) : std::vector<std::string>())
139 template<
typename EvalT,
typename Traits>
149 auto kokkosFieldMults_h = Kokkos::create_mirror_view(kokkosFieldMults_);
152 for (
size_t i(0); i < fieldMults_.size(); ++i)
153 kokkosFieldMults_h(i) = fieldMults_[i].get_static_view();
155 Kokkos::deep_copy(kokkosFieldMults_, kokkosFieldMults_h);
161 if (Sacado::IsADType<ScalarT>::value) {
162 const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
163 tmp_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",field_.extent(0),fadSize);
164 if (fieldMults_.size() > 1)
165 tmp2_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",field_.extent(0),fadSize);
167 tmp_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",field_.extent(0));
168 if (fieldMults_.size() > 1)
169 tmp2_ = PHX::View<ScalarT*>(
"IntegratorBasisTimesScalar::tmp_",field_.extent(0));
178 template<
typename EvalT,
typename Traits>
179 template<
int NUM_FIELD_MULT>
180 KOKKOS_INLINE_FUNCTION
185 const size_t& cell)
const
190 const int numQP(scalar_.extent(1)), numBases(basis_.extent(1));
192 for (
int basis(0); basis < numBases; ++basis)
193 field_(cell, basis) = 0.0;
197 if (NUM_FIELD_MULT == 0)
202 for (
int qp(0); qp < numQP; ++qp)
204 tmp_(cell) = multiplier_ * scalar_(cell, qp);
205 for (
int basis(0); basis < numBases; ++basis)
206 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
209 else if (NUM_FIELD_MULT == 1)
214 for (
int qp(0); qp < numQP; ++qp)
216 tmp_(cell) = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
217 for (
int basis(0); basis < numBases; ++basis)
218 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
227 const int numFieldMults(kokkosFieldMults_.extent(0));
228 for (
int qp(0); qp < numQP; ++qp)
231 for (
int fm(0); fm < numFieldMults; ++fm)
232 tmp2_(cell) *= kokkosFieldMults_(fm)(cell, qp);
233 tmp_(cell) = multiplier_ * scalar_(cell, qp) * tmp2_(cell);
234 for (
int basis(0); basis < numBases; ++basis)
235 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
245 template<
typename EvalT,
typename Traits>
251 using Kokkos::parallel_for;
252 using Kokkos::RangePolicy;
255 basis_ = this->wda(workset).bases[basisIndex_]->weighted_basis_scalar;
260 if (fieldMults_.size() == 0)
262 else if (fieldMults_.size() == 1)
273 template<
typename EvalT,
typename TRAITS>
288 p->set<
string>(
"Residual Name",
"?");
289 p->set<
string>(
"Value Name",
"?");
291 p->set(
"Basis", basis);
294 p->set<
double>(
"Multiplier", 1.0);
296 p->set(
"Field Multipliers", fms);
302 #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_