11 #ifndef __Panzer_Integerator_BasisTimesVector_impl_hpp__
12 #define __Panzer_Integerator_BasisTimesVector_impl_hpp__
33 template<
typename EvalT,
typename Traits>
37 const std::string& resName,
38 const std::string& valName,
42 const std::vector<std::string>& fmNames
45 evalStyle_(evalStyle),
46 useDescriptors_(false),
47 multiplier_(multiplier),
48 basisName_(basis.name())
58 using std::invalid_argument;
59 using std::logic_error;
65 "Integrator_BasisTimesVector called with an empty residual name.")
67 "Integrator_BasisTimesVector called with an empty value name.")
70 "Error: Integrator_BasisTimesVector: Basis of type \""
71 << tmpBasis->name() <<
"\" is not a vector basis.");
73 logic_error,
"Integrator_BasisTimesVector: Basis of type \""
74 << tmpBasis->name() <<
"\" does not require orientations. This seems " \
75 "very strange, so I'm failing.");
77 "ERROR: Integrator_BasisTimesVector only supports up to three multipliers!");
81 this->addDependentField(vector_);
87 this->addContributedField(
field_);
89 this->addEvaluatedField(
field_);
94 for (
const auto& name : fmNames) {
100 string n(
"Integrator_BasisTimesVector<");
101 n += std::to_string(fmNames.size()) +
">(";
106 n +=
", " + print<EvalT>() +
"): " +
field_.fieldTag().name();
115 template<
typename EvalT,
typename Traits>
124 const std::vector<PHX::FieldTag>& multipliers
127 evalStyle_(evalStyle),
128 useDescriptors_(true),
131 multiplier_(multiplier)
141 using std::invalid_argument;
142 using std::logic_error;
148 (
bd_.
getType() ==
"HDiv")), logic_error,
"Error: " \
149 "Integrator_BasisTimesVector: Basis of type \"" <<
bd_.
getType()
150 <<
"\" is not a vector basis.")
153 "ERROR: Integrator_BasisTimesVector only supports up to three multipliers!");
157 this->addDependentField(vector_);
163 this->addContributedField(
field_);
165 this->addEvaluatedField(
field_);
170 for (
const auto& fm : multipliers)
177 string n(
"Integrator_BasisTimesVector<");
178 n += std::to_string(multipliers.size()) +
">(";
183 n +=
", " + print<EvalT>() +
"): " +
field_.fieldTag().name();
192 template<
typename EvalT,
typename Traits>
199 p.get<std::string>(
"Residual Name"),
200 p.get<std::string>(
"Value Name"),
203 p.get<double>(
"Multiplier"),
204 p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
205 (
"Field Multipliers") ?
206 (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
207 (
"Field Multipliers")) : std::vector<std::string>())
222 template<
typename EvalT,
typename Traits>
233 for (
size_t i=0; i < fieldMults_.size(); ++i)
234 kokkosFieldMults_[i] = fieldMults_[i].get_static_view();
238 numQP_ = vector_.extent(1);
239 numDim_ = vector_.extent(2);
242 if (not useDescriptors_)
252 template<
typename EvalT,
typename Traits>
253 template<
int NUM_FIELD_MULT>
254 KOKKOS_INLINE_FUNCTION
262 const int cell = team.league_rank();
263 const int numBases = basis_.extent(1);
267 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
268 field_(cell, basis) = 0.0;
272 for (
int qp(0); qp < numQP_; ++qp) {
273 for (
int dim(0); dim < numDim_; ++dim) {
274 if constexpr (NUM_FIELD_MULT == 0) {
275 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
276 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim);
279 else if constexpr (NUM_FIELD_MULT == 1) {
280 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
281 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
282 * kokkosFieldMults_[0](cell, qp);
285 else if constexpr (NUM_FIELD_MULT == 2) {
286 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
287 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
288 * kokkosFieldMults_[0](cell, qp) * kokkosFieldMults_[1](cell, qp);
291 else if constexpr (NUM_FIELD_MULT == 3) {
292 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
293 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
294 * kokkosFieldMults_[0](cell, qp) * kokkosFieldMults_[1](cell, qp) * kokkosFieldMults_[2](cell, qp);
298 Kokkos::abort(
"Panzer_Integrator_BasisTimesVector: NUM_FIELD_MULT out of range!");
309 template<
typename EvalT,
typename Traits>
315 using Kokkos::parallel_for;
316 using Kokkos::RangePolicy;
321 *this->wda(workset).bases[basisIndex_];
328 if (fieldMults_.size() == 0) {
330 parallel_for(
"Panzer_Integrator_BasisTimesVector<0>", policy, *
this);
332 else if (fieldMults_.size() == 1) {
334 parallel_for(
"Panzer_Integrator_BasisTimesVector<1>", policy, *
this);
336 else if (fieldMults_.size() == 2) {
338 parallel_for(
"Panzer_Integrator_BasisTimesVector<2>", policy, *
this);
340 else if (fieldMults_.size() == 3) {
342 parallel_for(
"Panzer_Integrator_BasisTimesVector<3>", policy, *
this);
351 template<
typename EvalT,
typename Traits>
366 p->set<
string>(
"Residual Name",
"?");
367 p->set<
string>(
"Value Name",
"?");
369 p->set(
"Basis", basis);
372 p->set<
double>(
"Multiplier", 1.0);
374 p->set(
"Field Multipliers", fms);
380 #endif // __Panzer_Integerator_BasisTimesVector_impl_hpp__
int num_cells
DEPRECATED - use: numCells()
Kokkos::TeamPolicy< TeamPolicyProperties...> teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
BasisDescriptor bd_
The BasisDescriptor for the basis to use.
EvaluatorStyle
An indication of how an Evaluator will behave.
Teuchos::RCP< const PureBasis > getBasis() const
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
double multiplier
The scalar multiplier out in front of the integral ( ).
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const typename Kokkos::TeamPolicy< FieldMultTag< NUM_FIELD_MULT >, PHX::exec_space >::member_type &team) const
Perform the integration.
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_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const std::string & getType() const
Get type of basis.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
typename EvalT::ScalarT ScalarT
The scalar type.
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.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
Array_CellBasisIPDim weighted_basis_vector
ConstArray_CellBasisIPDim getVectorBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at mesh points.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
static HP & inst()
Private ctor.
Integrator_BasisTimesVector(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.
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > vector_
A field representing the vector-valued function we're integrating ( ).
Description and data layouts associated with a particular 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< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_