43 #ifndef __Panzer_Integrator_GradBasisDotVector_impl_hpp__
44 #define __Panzer_Integrator_GradBasisDotVector_impl_hpp__
65 template<
typename EvalT,
typename Traits>
69 const std::string& resName,
70 const std::string& fluxName,
74 const std::vector<std::string>& fmNames,
78 evalStyle_(evalStyle),
79 multiplier_(multiplier),
80 basisName_(basis.name())
89 using std::invalid_argument;
90 using std::logic_error;
96 "Integrator_GradBasisDotVector called with an empty residual name.")
98 "Integrator_GradBasisDotVector called with an empty flux name.")
101 "Error: Integrator_GradBasisDotVector: Basis of type \""
102 << tmpBasis->name() <<
"\" does not support the gradient operator.")
108 tmpVecDL->extent(2) < ir.
dl_vector->extent(2), logic_error,
109 "Integrator_GradBasisDotVector: Dimension of space exceeds " \
110 "dimension of Vector Data Layout.");
114 vector_ = MDField<const ScalarT, Cell, IP, Dim>(fluxName, tmpVecDL);
115 this->addDependentField(
vector_);
121 this->addContributedField(
field_);
123 this->addEvaluatedField(
field_);
128 kokkosFieldMults_ = PHX::View<PHX::UnmanagedView<const ScalarT**>*>(
"GradBasisDotVector::KokkosFieldMultipliers",
130 for (
const auto& name : fmNames)
137 string n(
"Integrator_GradBasisDotVector (");
142 n +=
"): " +
field_.fieldTag().name();
151 template<
typename EvalT,
typename Traits>
158 p.get<std::string>(
"Residual Name"),
159 p.get<std::string>(
"Flux Name"),
162 p.get<double>(
"Multiplier"),
163 p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
164 (
"Field Multipliers") ?
165 (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
166 (
"Field Multipliers")) : std::vector<std::string>(),
167 p.isType<Teuchos::
RCP<PHX::DataLayout>>(
"Vector Data Layout") ?
168 p.get<Teuchos::
RCP<PHX::DataLayout>>(
"Vector Data Layout") :
184 template<
typename EvalT,
typename Traits>
195 auto field_mults_host_mirror = Kokkos::create_mirror_view(kokkosFieldMults_);
196 for (
size_t i(0); i < fieldMults_.size(); ++i)
197 field_mults_host_mirror(i) = fieldMults_[i].get_static_view();
198 Kokkos::deep_copy(kokkosFieldMults_,field_mults_host_mirror);
205 if (!use_shared_memory) {
206 if (Sacado::IsADType<ScalarT>::value) {
207 const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
208 tmp_ = PHX::View<ScalarT*>(
"GradBasisDotVector::tmp_",field_.extent(0),fadSize);
210 tmp_ = PHX::View<ScalarT*>(
"GradBasisDotVector::tmp_",field_.extent(0));
220 template<
typename EvalT,
typename Traits>
221 template<
int NUM_FIELD_MULT>
222 KOKKOS_INLINE_FUNCTION
227 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
230 const int cell = team.league_rank();
233 const int numQP(vector_.extent(1)), numDim(vector_.extent(2)),
234 numBases(basis_.extent(1));
236 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
237 field_(cell, basis) = 0.0;
242 if (NUM_FIELD_MULT == 0)
247 for (
int qp(0); qp < numQP; ++qp)
249 for (
int dim(0); dim < numDim; ++dim)
251 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
252 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim);
257 else if (NUM_FIELD_MULT == 1)
262 for (
int qp(0); qp < numQP; ++qp)
264 for (
int dim(0); dim < numDim; ++dim)
266 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
267 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim) * kokkosFieldMults_(0)(cell, qp);
279 const int numFieldMults(kokkosFieldMults_.extent(0));
280 for (
int qp(0); qp < numQP; ++qp)
284 for (
int fm(0); fm < numFieldMults; ++fm)
285 tmp_(cell) *= kokkosFieldMults_(fm)(cell, qp);
286 for (
int dim(0); dim < numDim; ++dim)
288 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
289 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim) * tmp_(cell);
301 template<
typename EvalT,
typename Traits>
302 template<
int NUM_FIELD_MULT>
303 KOKKOS_INLINE_FUNCTION
308 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
311 const int cell = team.league_rank();
312 const int numQP = vector_.extent(1);
313 const int numDim = vector_.extent(2);
314 const int numBases = basis_.extent(1);
315 const int fadSize = Kokkos::dimension_scalar(field_.get_view());
319 if (Sacado::IsADType<ScalarT>::value) {
321 tmp_field =
scratch_view(team.team_shmem(),numBases,fadSize);
329 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
330 tmp_field(basis) = 0.0;
335 if (NUM_FIELD_MULT == 0)
340 for (
int qp(0); qp < numQP; ++qp)
342 for (
int dim(0); dim < numDim; ++dim)
344 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
345 tmp_field(basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim);
350 else if (NUM_FIELD_MULT == 1)
355 for (
int qp(0); qp < numQP; ++qp)
357 for (
int dim(0); dim < numDim; ++dim)
359 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
360 tmp_field(basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim) * kokkosFieldMults_(0)(cell, qp);
372 const int numFieldMults(kokkosFieldMults_.extent(0));
373 for (
int qp(0); qp < numQP; ++qp)
376 for (
int fm(0); fm < numFieldMults; ++fm)
377 fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
378 for (
int dim(0); dim < numDim; ++dim)
381 tmp(0) = multiplier_ * vector_(cell, qp, dim) * fieldMultsTotal;
382 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
383 tmp_field(basis) += basis_(cell, basis, qp, dim) * tmp(0);
391 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int basis) {
392 field_(cell,basis) = tmp_field(basis);
396 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int basis) {
397 field_(cell,basis) += tmp_field(basis);
408 template<
typename EvalT,
typename Traits>
414 using Kokkos::parallel_for;
415 using Kokkos::TeamPolicy;
418 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
421 if (use_shared_memory) {
423 if (Sacado::IsADType<ScalarT>::value) {
424 const int fadSize = Kokkos::dimension_scalar(field_.get_view());
425 bytes = scratch_view::shmem_size(1,fadSize) + scratch_view::shmem_size(basis_.extent(1),fadSize);
428 bytes = scratch_view::shmem_size(1) + scratch_view::shmem_size(basis_.extent(1));
433 if (fieldMults_.size() == 0) {
435 parallel_for(this->getName(), policy, *
this);
436 }
else if (fieldMults_.size() == 1) {
438 parallel_for(this->getName(), policy, *
this);
441 parallel_for(this->getName(), policy, *
this);
448 if (fieldMults_.size() == 0) {
450 parallel_for(this->getName(), policy, *
this);
451 }
else if (fieldMults_.size() == 1) {
453 parallel_for(this->getName(), policy, *
this);
456 parallel_for(this->getName(), policy, *
this);
466 template<
typename EvalT,
typename TRAITS>
482 p->set<
string>(
"Residual Name",
"?");
483 p->set<
string>(
"Flux Name",
"?");
485 p->set(
"Basis", basis);
488 p->set<
double>(
"Multiplier", 1.0);
490 p->set(
"Field Multipliers", fms);
492 p->set(
"Vector Data Layout", vecDL);
498 #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_