43 #ifndef __Panzer_Integrator_DivBasisTimesScalar_impl_hpp__ 44 #define __Panzer_Integrator_DivBasisTimesScalar_impl_hpp__ 53 #include "Intrepid2_FunctionSpaceTools.hpp" 56 #include "Kokkos_ViewFactory.hpp" 71 template<
typename EvalT,
typename Traits>
75 const std::string& resName,
76 const std::string& valName,
80 const std::vector<std::string>& fmNames
85 basisName_(basis.name()),
86 use_shared_memory(false)
95 using std::invalid_argument;
96 using std::logic_error;
101 TEUCHOS_TEST_FOR_EXCEPTION(resName ==
"", invalid_argument,
"Error: " \
102 "Integrator_DivBasisTimesScalar called with an empty residual name.")
103 TEUCHOS_TEST_FOR_EXCEPTION(valName ==
"", invalid_argument,
"Error: " \
104 "Integrator_DivBasisTimesScalar called with an empty value name.")
105 RCP<const PureBasis> tmpBasis = basis.
getBasis();
106 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsDiv(), logic_error,
107 "Error: Integrator_DivBasisTimesScalar: Basis of type \"" 108 << tmpBasis->name() <<
"\" does not support DIV.")
109 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(), logic_error,
110 "Error: Integration_DivBasisTimesScalar: Basis of type \"" 111 << tmpBasis->name() <<
"\" should require orientations.")
115 this->addDependentField(
scalar_);
121 this->addContributedField(
field_);
123 this->addEvaluatedField(
field_);
129 View<View<const ScalarT**>*>(
"BasisTimesScalar::KokkosFieldMultipliers",
131 for (
const auto& name : fmNames)
138 string n(
"Integrator_DivBasisTimesScalar (");
143 n +=
"): " +
field_.fieldTag().name();
152 template<
typename EvalT,
typename Traits>
155 const Teuchos::ParameterList& p)
159 p.get<
std::string>(
"Residual Name"),
160 p.get<
std::string>(
"Value Name"),
163 p.get<double>(
"Multiplier"),
165 (
"Field Multipliers") ?
167 (
"Field Multipliers")) :
std::vector<
std::string>())
169 using Teuchos::ParameterList;
174 p.validateParameters(*validParams);
182 template<
typename EvalT,
typename Traits>
193 for (
size_t i(0); i < fieldMults_.size(); ++i)
194 kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
205 template<
typename EvalT,
typename Traits>
206 template<
int NUM_FIELD_MULT>
207 KOKKOS_INLINE_FUNCTION
212 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const 215 const int cell = team.league_rank();
218 const int numQP(scalar_.extent(1)), numBases(basis_.extent(1));
220 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
221 field_(cell, basis) = 0.0;
228 if (NUM_FIELD_MULT == 0)
233 for (
int qp(0); qp < numQP; ++qp)
235 tmp = multiplier_ * scalar_(cell, qp);
236 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
237 field_(cell, basis) += basis_(cell, basis, qp) * tmp;
241 else if (NUM_FIELD_MULT == 1)
246 for (
int qp(0); qp < numQP; ++qp)
248 tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
249 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
250 field_(cell, basis) += basis_(cell, basis, qp) * tmp;
260 const int numFieldMults(kokkosFieldMults_.extent(0));
261 for (
int qp(0); qp < numQP; ++qp)
264 for (
int fm(0); fm < numFieldMults; ++fm)
265 fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
266 tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
267 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
268 field_(cell, basis) += basis_(cell, basis, qp) * tmp;
279 template<
typename EvalT,
typename Traits>
280 template<
int NUM_FIELD_MULT>
281 KOKKOS_INLINE_FUNCTION
286 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const 289 const int cell = team.league_rank();
290 const int numQP = scalar_.extent(1);
291 const int numBases = basis_.extent(1);
292 const int fadSize = Kokkos::dimension_scalar(field_.get_view());
296 if (Sacado::IsADType<ScalarT>::value) {
298 tmp_field =
scratch_view(team.team_shmem(),numBases,fadSize);
306 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
307 tmp_field(basis) = 0.0;
312 if (NUM_FIELD_MULT == 0)
317 for (
int qp(0); qp < numQP; ++qp)
320 tmp(0) = multiplier_ * scalar_(cell, qp);
321 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
322 tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
326 else if (NUM_FIELD_MULT == 1)
331 for (
int qp(0); qp < numQP; ++qp)
334 tmp(0) = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
335 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
336 tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
346 const int numFieldMults = kokkosFieldMults_.extent(0);
347 for (
int qp(0); qp < numQP; ++qp)
351 for (
int fm(0); fm < numFieldMults; ++fm)
352 fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
353 tmp(0) = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
354 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
355 tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
362 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int& basis) {
363 field_(cell,basis) = tmp_field(basis);
367 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int& basis) {
368 field_(cell,basis) += tmp_field(basis);
379 template<
typename EvalT,
typename Traits>
385 using Kokkos::parallel_for;
386 using Kokkos::TeamPolicy;
389 basis_ = this->wda(workset).bases[basisIndex_]->weighted_div_basis;
393 if (use_shared_memory) {
395 if (Sacado::IsADType<ScalarT>::value) {
396 const int fadSize = Kokkos::dimension_scalar(field_.get_view());
397 bytes = scratch_view::shmem_size(1,fadSize) + scratch_view::shmem_size(basis_.extent(1),fadSize);
400 bytes = scratch_view::shmem_size(1) + scratch_view::shmem_size(basis_.extent(1));
405 if (fieldMults_.size() == 0) {
407 parallel_for(policy, *
this, this->getName());
408 }
else if (fieldMults_.size() == 1) {
410 parallel_for(policy, *
this, this->getName());
413 parallel_for(policy, *
this, this->getName());
420 if (fieldMults_.size() == 0) {
422 parallel_for(policy, *
this, this->getName());
423 }
else if (fieldMults_.size() == 1) {
425 parallel_for(policy, *
this, this->getName());
428 parallel_for(policy, *
this, this->getName());
438 template<
typename EvalT,
typename TRAITS>
439 Teuchos::RCP<Teuchos::ParameterList>
446 using Teuchos::ParameterList;
451 RCP<ParameterList> p = rcp(
new ParameterList);
452 p->set<
string>(
"Residual Name",
"?");
453 p->set<
string>(
"Value Name",
"?");
454 p->set<
string>(
"Test Field Name",
"?");
455 RCP<BasisIRLayout> basis;
456 p->set(
"Basis", basis);
457 RCP<IntegrationRule> ir;
459 p->set<
double>(
"Multiplier", 1.0);
460 RCP<const vector<string>> fms;
461 p->set(
"Field Multipliers", fms);
467 #endif // __Panzer_Integrator_DivBasisTimesScalar_impl_hpp__ int num_cells
DEPRECATED - use: numCells()
bool useSharedMemory() const
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
Kokkos::TeamPolicy< TeamPolicyProperties... > teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
typename EvalT::ScalarT ScalarT
The scalar type.
PHX::View< PHX::View< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
Teuchos::RCP< const PureBasis > getBasis() const
EvaluatorStyle
An indication of how an Evaluator will behave.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
double multiplier
The scalar multiplier out in front of the integral ( ).
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we're integrating ( ).
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack... dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
Perform the integration. Main memory version.
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< ScalarT, Cell, BASIS > field_
A field to which we'll contribute, or in which we'll store, the result of computing this integral...
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Integrator_DivBasisTimesScalar(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.
static HP & inst()
Private ctor.
Kokkos::View< ScalarT *,typename PHX::DevLayout< ScalarT >::type, typename PHX::exec_space::scratch_memory_space, Kokkos::MemoryUnmanaged > scratch_view
Type for shared memory.
Description and data layouts associated with a particular basis.
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_
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ...