11 #ifndef __Panzer_Integrator_GradBasisDotVector_impl_hpp__
12 #define __Panzer_Integrator_GradBasisDotVector_impl_hpp__
33 template<
typename EvalT,
typename Traits>
37 const std::string& resName,
38 const std::string& fluxName,
42 const std::vector<std::string>& fmNames,
46 evalStyle_(evalStyle),
47 multiplier_(multiplier),
48 basisName_(basis.name())
57 using std::invalid_argument;
58 using std::logic_error;
64 "Integrator_GradBasisDotVector called with an empty residual name.")
66 "Integrator_GradBasisDotVector called with an empty flux name.")
69 "Error: Integrator_GradBasisDotVector: Basis of type \""
70 << tmpBasis->name() <<
"\" does not support the gradient operator.")
76 tmpVecDL->extent(2) < ir.
dl_vector->extent(2), logic_error,
77 "Integrator_GradBasisDotVector: Dimension of space exceeds " \
78 "dimension of Vector Data Layout.");
82 vector_ = MDField<const ScalarT, Cell, IP, Dim>(fluxName, tmpVecDL);
83 this->addDependentField(
vector_);
89 this->addContributedField(
field_);
91 this->addEvaluatedField(
field_);
96 kokkosFieldMults_ = PHX::View<PHX::UnmanagedView<const ScalarT**>*>(
"GradBasisDotVector::KokkosFieldMultipliers",
98 for (
const auto& name : fmNames)
105 string n(
"Integrator_GradBasisDotVector (");
110 n +=
"): " +
field_.fieldTag().name();
119 template<
typename EvalT,
typename Traits>
126 p.get<std::string>(
"Residual Name"),
127 p.get<std::string>(
"Flux Name"),
130 p.get<double>(
"Multiplier"),
131 p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
132 (
"Field Multipliers") ?
133 (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
134 (
"Field Multipliers")) : std::vector<std::string>(),
135 p.isType<Teuchos::
RCP<PHX::DataLayout>>(
"Vector Data Layout") ?
136 p.get<Teuchos::
RCP<PHX::DataLayout>>(
"Vector Data Layout") :
152 template<
typename EvalT,
typename Traits>
163 auto field_mults_host_mirror = Kokkos::create_mirror_view(kokkosFieldMults_);
164 for (
size_t i(0); i < fieldMults_.size(); ++i)
165 field_mults_host_mirror(i) = fieldMults_[i].get_static_view();
166 Kokkos::deep_copy(kokkosFieldMults_,field_mults_host_mirror);
173 if (!use_shared_memory) {
174 if (Sacado::IsADType<ScalarT>::value) {
175 const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
176 tmp_ = PHX::View<ScalarT*>(
"GradBasisDotVector::tmp_",field_.extent(0),fadSize);
178 tmp_ = PHX::View<ScalarT*>(
"GradBasisDotVector::tmp_",field_.extent(0));
188 template<
typename EvalT,
typename Traits>
189 template<
int NUM_FIELD_MULT>
190 KOKKOS_INLINE_FUNCTION
195 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
198 const int cell = team.league_rank();
201 const int numQP(vector_.extent(1)), numDim(vector_.extent(2)),
202 numBases(basis_.extent(1));
204 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
205 field_(cell, basis) = 0.0;
208 for (
int qp(0); qp < numQP; ++qp) {
209 for (
int dim(0); dim < numDim; ++dim) {
210 if constexpr (NUM_FIELD_MULT == 0) {
211 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
212 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim);
215 else if constexpr (NUM_FIELD_MULT == 1) {
216 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
217 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
218 * kokkosFieldMults_(0)(cell, qp);
221 else if constexpr (NUM_FIELD_MULT == 2) {
222 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
223 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
224 * kokkosFieldMults_(0)(cell, qp) * kokkosFieldMults_(1)(cell, qp);
227 else if constexpr (NUM_FIELD_MULT == 3) {
228 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
229 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
230 * kokkosFieldMults_(0)(cell, qp) * kokkosFieldMults_(1)(cell, qp) * kokkosFieldMults_(2)(cell, qp);
234 Kokkos::abort(
"Panzer_Integrator_GradBasisDotVector: NUM_FIELD_MULT out of range!");
245 template<
typename EvalT,
typename Traits>
246 template<
int NUM_FIELD_MULT>
247 KOKKOS_INLINE_FUNCTION
252 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
255 const int cell = team.league_rank();
256 const int numQP = vector_.extent(1);
257 const int numDim = vector_.extent(2);
258 const int numBases = basis_.extent(1);
259 const int fadSize = Kokkos::dimension_scalar(field_.get_view());
262 if (Sacado::IsADType<ScalarT>::value) {
263 tmp_field =
scratch_view(team.team_shmem(),numBases,fadSize);
270 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
271 tmp_field(basis) = 0.0;
276 for (
int qp(0); qp < numQP; ++qp) {
277 for (
int dim(0); dim < numDim; ++dim) {
278 if constexpr (NUM_FIELD_MULT == 0) {
279 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
280 tmp_field(basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim);
283 else if constexpr (NUM_FIELD_MULT == 1) {
284 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
285 tmp_field(basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
286 * kokkosFieldMults_(0)(cell, qp);
289 else if constexpr (NUM_FIELD_MULT == 2) {
290 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
291 tmp_field(basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
292 * kokkosFieldMults_(0)(cell, qp) * kokkosFieldMults_(1)(cell, qp);
295 else if constexpr (NUM_FIELD_MULT == 3) {
296 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
297 tmp_field(basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
298 * kokkosFieldMults_(0)(cell, qp) * kokkosFieldMults_(1)(cell, qp) * kokkosFieldMults_(2)(cell, qp);
302 Kokkos::abort(
"Panzer_Integrator_GradBasisDotVector: NUM_FIELD_MULT out of range!");
309 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int basis) {
310 field_(cell,basis) = tmp_field(basis);
314 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int basis) {
315 field_(cell,basis) += tmp_field(basis);
326 template<
typename EvalT,
typename Traits>
332 using Kokkos::parallel_for;
333 using Kokkos::TeamPolicy;
336 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
339 if (use_shared_memory) {
341 if (Sacado::IsADType<ScalarT>::value) {
342 const int fadSize = Kokkos::dimension_scalar(field_.get_view());
343 bytes = scratch_view::shmem_size(basis_.extent(1),fadSize);
346 bytes = scratch_view::shmem_size(basis_.extent(1));
351 if (fieldMults_.size() == 0) {
353 parallel_for(this->getName(), policy, *
this);
354 }
else if (fieldMults_.size() == 1) {
356 parallel_for(this->getName(), policy, *
this);
357 }
else if (fieldMults_.size() == 2) {
359 parallel_for(this->getName(), policy, *
this);
360 }
else if (fieldMults_.size() == 3) {
362 parallel_for(this->getName(), policy, *
this);
365 "ERROR: Panzer_Integrator_GradBasisDotVector supports up to three field multipliers! User requested "
366 << fieldMults_.size() <<
"!");
373 if (fieldMults_.size() == 0) {
375 parallel_for(this->getName(), policy, *
this);
376 }
else if (fieldMults_.size() == 1) {
378 parallel_for(this->getName(), policy, *
this);
379 }
else if (fieldMults_.size() == 2) {
381 parallel_for(this->getName(), policy, *
this);
382 }
else if (fieldMults_.size() == 3) {
384 parallel_for(this->getName(), policy, *
this);
387 "ERROR: Panzer_Integrator_GradBasisDotVector supports up to three field multipliers! User requested "
388 << fieldMults_.size() <<
"!");
398 template<
typename EvalT,
typename TRAITS>
414 p->set<
string>(
"Residual Name",
"?");
415 p->set<
string>(
"Flux Name",
"?");
417 p->set(
"Basis", basis);
420 p->set<
double>(
"Multiplier", 1.0);
422 p->set(
"Field Multipliers", fms);
424 p->set(
"Vector Data Layout", vecDL);
430 #endif // __Panzer_Integrator_GradBasisDotVector_impl_hpp__
int num_cells
DEPRECATED - use: numCells()
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > vector_
A field representing the vector-valued function we're integrating ( ).
typename EvalT::ScalarT ScalarT
The scalar type.
Integrator_GradBasisDotVector(const panzer::EvaluatorStyle &evalStyle, const std::string &resName, const std::string &fluxName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >(), const Teuchos::RCP< PHX::DataLayout > &vecDL=Teuchos::null)
Main Constructor.
Kokkos::TeamPolicy< TeamPolicyProperties...> teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
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)
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
Perform the integration.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
EvaluatorStyle
An indication of how an Evaluator will behave.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
Teuchos::RCP< const PureBasis > getBasis() const
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 ( ...
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.
bool useSharedMemory() const
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.
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
static HP & inst()
Private ctor.
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_