43 #ifndef __Panzer_Integerator_BasisTimesVector_impl_hpp__
44 #define __Panzer_Integerator_BasisTimesVector_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 useDescriptors_(false),
78 multiplier_(multiplier),
79 basisName_(basis.name())
89 using std::invalid_argument;
90 using std::logic_error;
96 "Integrator_BasisTimesVector called with an empty residual name.")
98 "Integrator_BasisTimesVector called with an empty value name.")
101 "Error: Integrator_BasisTimesVector: Basis of type \""
102 << tmpBasis->name() <<
"\" is not a vector basis.");
104 logic_error,
"Integrator_BasisTimesVector: Basis of type \""
105 << tmpBasis->name() <<
"\" does not require orientations. This seems " \
106 "very strange, so I'm failing.");
110 this->addDependentField(vector_);
116 this->addContributedField(
field_);
118 this->addEvaluatedField(
field_);
124 View<View<const ScalarT**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>(
"BasisTimesVector::KokkosFieldMultipliers",
126 for (
const auto& name : fmNames)
133 string n(
"Integrator_BasisTimesVector (");
138 n +=
", " + print<EvalT>() +
"): " +
field_.fieldTag().name();
147 template<
typename EvalT,
typename Traits>
151 const PHX::FieldTag& resTag,
152 const PHX::FieldTag& valTag,
156 const std::vector<PHX::FieldTag>& multipliers
159 evalStyle_(evalStyle),
160 useDescriptors_(true),
163 multiplier_(multiplier)
173 using std::invalid_argument;
174 using std::logic_error;
180 (
bd_.
getType() ==
"HDiv")), logic_error,
"Error: " \
181 "Integrator_BasisTimesVector: Basis of type \"" <<
bd_.
getType()
182 <<
"\" is not a vector basis.")
186 this->addDependentField(
vector_);
192 this->addContributedField(
field_);
194 this->addEvaluatedField(
field_);
200 View<View<const ScalarT**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>(
"BasisTimesVector::KokkosFieldMultipliers",
202 for (
const auto& fm : multipliers)
209 string n(
"Integrator_BasisTimesVector (");
214 n +=
", " + print<EvalT>() +
"): " +
field_.fieldTag().name();
223 template<
typename EvalT,
typename Traits>
230 p.get<std::string>(
"Residual Name"),
231 p.get<std::string>(
"Value Name"),
234 p.get<double>(
"Multiplier"),
235 p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
236 (
"Field Multipliers") ?
237 (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
238 (
"Field Multipliers")) : std::vector<std::string>())
253 template<
typename EvalT,
typename Traits>
264 for (
size_t i(0); i < fieldMults_.size(); ++i)
265 kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
269 numQP_ = vector_.extent(1);
270 numDim_ = vector_.extent(2);
273 if (not useDescriptors_)
282 template<
typename EvalT,
typename Traits>
283 template<
int NUM_FIELD_MULT>
284 KOKKOS_INLINE_FUNCTION
289 const size_t& cell)
const
294 const int numBases(basis_.extent(1));
296 for (
int basis(0); basis < numBases; ++basis)
297 field_(cell, basis) = 0.0;
302 if (NUM_FIELD_MULT == 0)
307 for (
int qp(0); qp < numQP_; ++qp)
309 for (
int dim(0); dim < numDim_; ++dim)
311 tmp = multiplier_ * vector_(cell, qp, dim);
312 for (
int basis(0); basis < numBases; ++basis)
313 field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
317 else if (NUM_FIELD_MULT == 1)
322 for (
int qp(0); qp < numQP_; ++qp)
324 for (
int dim(0); dim < numDim_; ++dim)
326 tmp = multiplier_ * vector_(cell, qp, dim) *
327 kokkosFieldMults_(0)(cell, qp);
328 for (
int basis(0); basis < numBases; ++basis)
329 field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
340 const int numFieldMults(kokkosFieldMults_.extent(0));
341 for (
int qp(0); qp < numQP_; ++qp)
344 for (
int fm(0); fm < numFieldMults; ++fm)
345 fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
346 for (
int dim(0); dim < numDim_; ++dim)
348 tmp = multiplier_ * vector_(cell, qp, dim) * fieldMultsTotal;
349 for (
int basis(0); basis < numBases; ++basis)
350 field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
361 template<
typename EvalT,
typename Traits>
367 using Kokkos::parallel_for;
368 using Kokkos::RangePolicy;
372 this->wda(workset).getBasisValues(bd_,id_) :
373 *this->wda(workset).bases[basisIndex_];
379 if (fieldMults_.size() == 0)
381 else if (fieldMults_.size() == 1)
392 template<
typename EvalT,
typename Traits>
407 p->set<
string>(
"Residual Name",
"?");
408 p->set<
string>(
"Value Name",
"?");
410 p->set(
"Basis", basis);
413 p->set<
double>(
"Multiplier", 1.0);
415 p->set(
"Field Multipliers", fms);
421 #endif // __Panzer_Integerator_BasisTimesVector_impl_hpp__
#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 ( ).
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.
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
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
Kokkos::View< Kokkos::View< const ScalarT **, typename PHX::DevLayout< ScalarT >::type, PHX::Device > * > kokkosFieldMults_
The Kokkos::View representation of the (possibly empty) list of fields that are multipliers out in fr...
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>
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_