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     
 
   83     evalStyle_(evalStyle),
 
   84     multiplier_(multiplier),
 
   85     basisName_(basis.name()),
 
   86     use_shared_memory(false)
 
   95     using std::invalid_argument;
 
   96     using std::logic_error;
 
  102       "Integrator_DivBasisTimesScalar called with an empty residual name.")
 
  104       "Integrator_DivBasisTimesScalar called with an empty value name.")
 
  107       "Error:  Integrator_DivBasisTimesScalar:  Basis of type \"" 
  108       << tmpBasis->name() << 
"\" does not support DIV.")
 
  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**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>(
"BasisTimesScalar::KokkosFieldMultipliers",
 
  131     for (
const auto& name : fmNames)
 
  138     string n(
"Integrator_DivBasisTimesScalar (");
 
  143     n += 
"):  " + 
field_.fieldTag().name();
 
  152   template<
typename EvalT, 
typename Traits>
 
  159       p.get<std::string>(
"Residual Name"),
 
  160       p.get<std::string>(
"Value Name"),
 
  163       p.get<double>(
"Multiplier"),
 
  164       p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
 
  165         (
"Field Multipliers") ?
 
  166         (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
 
  167         (
"Field Multipliers")) : std::vector<std::string>())
 
  182   template<
typename EvalT, 
typename Traits>
 
  194     for (
size_t i(0); i < fieldMults_.size(); ++i)
 
  195       kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
 
  207   template<
typename EvalT, 
typename Traits>
 
  208   template<
int NUM_FIELD_MULT>
 
  209   KOKKOS_INLINE_FUNCTION
 
  214     const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
 const 
  217     const int cell = team.league_rank();
 
  220     const int numQP(scalar_.extent(1)), numBases(basis_.extent(1));
 
  222       Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  223   field_(cell, basis) = 0.0;
 
  230     if (NUM_FIELD_MULT == 0)
 
  235       for (
int qp(0); qp < numQP; ++qp)
 
  237         tmp = multiplier_ * scalar_(cell, qp);
 
  238   Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  239           field_(cell, basis) += basis_(cell, basis, qp) * tmp;
 
  243     else if (NUM_FIELD_MULT == 1)
 
  248       for (
int qp(0); qp < numQP; ++qp)
 
  250         tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
 
  251   Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  252           field_(cell, basis) += basis_(cell, basis, qp) * tmp;
 
  262       const int numFieldMults(kokkosFieldMults_.extent(0));
 
  263       for (
int qp(0); qp < numQP; ++qp)
 
  266         for (
int fm(0); fm < numFieldMults; ++fm)
 
  267           fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
 
  268         tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
 
  269   Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  270           field_(cell, basis) += basis_(cell, basis, qp) * tmp;
 
  281   template<
typename EvalT, 
typename Traits>
 
  282   template<
int NUM_FIELD_MULT>
 
  283   KOKKOS_INLINE_FUNCTION
 
  288     const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
 const 
  291     const int cell = team.league_rank();
 
  292     const int numQP = scalar_.extent(1);
 
  293     const int numBases = basis_.extent(1);
 
  294     const int fadSize = Kokkos::dimension_scalar(field_.get_view());
 
  298     if (Sacado::IsADType<ScalarT>::value) {
 
  300       tmp_field = 
scratch_view(team.team_shmem(),numBases,fadSize);
 
  308     Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  309       tmp_field(basis) = 0.0;
 
  314     if (NUM_FIELD_MULT == 0)
 
  319       for (
int qp(0); qp < numQP; ++qp)
 
  322   tmp(0) = multiplier_ * scalar_(cell, qp);
 
  323   Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  324     tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
 
  328     else if (NUM_FIELD_MULT == 1)
 
  333       for (
int qp(0); qp < numQP; ++qp)
 
  336         tmp(0) = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
 
  337   Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  338     tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
 
  348       const int numFieldMults = kokkosFieldMults_.extent(0);
 
  349       for (
int qp(0); qp < numQP; ++qp)
 
  353         for (
int fm(0); fm < numFieldMults; ++fm)
 
  354           fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
 
  355         tmp(0) = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
 
  356   Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  357     tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
 
  364       Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int& basis) {
 
  365   field_(cell,basis) = tmp_field(basis);
 
  369       Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int& basis) {
 
  370   field_(cell,basis) += tmp_field(basis);
 
  381   template<
typename EvalT, 
typename Traits>
 
  387     using Kokkos::parallel_for;
 
  388     using Kokkos::TeamPolicy;
 
  391     basis_ = this->wda(workset).bases[basisIndex_]->weighted_div_basis;
 
  395     if (use_shared_memory) {
 
  397       if (Sacado::IsADType<ScalarT>::value) {
 
  398   const int fadSize = Kokkos::dimension_scalar(field_.get_view());
 
  399   bytes = scratch_view::shmem_size(1,fadSize) + scratch_view::shmem_size(basis_.extent(1),fadSize);
 
  402   bytes = scratch_view::shmem_size(1) + scratch_view::shmem_size(basis_.extent(1));
 
  407       if (fieldMults_.size() == 0) {
 
  409   parallel_for(policy, *
this, this->getName());
 
  410       } 
else if (fieldMults_.size() == 1) {
 
  412   parallel_for(policy, *
this, this->getName());
 
  415   parallel_for(policy, *
this, this->getName());
 
  422       if (fieldMults_.size() == 0) {
 
  424   parallel_for(policy, *
this, this->getName());
 
  425       } 
else if (fieldMults_.size() == 1) {
 
  427   parallel_for(policy, *
this, this->getName());
 
  430   parallel_for(policy, *
this, this->getName());
 
  440   template<
typename EvalT, 
typename TRAITS>
 
  454     p->set<
string>(
"Residual Name", 
"?");
 
  455     p->set<
string>(
"Value Name", 
"?");
 
  456     p->set<
string>(
"Test Field Name", 
"?");
 
  458     p->set(
"Basis", basis);
 
  461     p->set<
double>(
"Multiplier", 1.0);
 
  463     p->set(
"Field Multipliers", fms);
 
  469 #endif // __Panzer_Integrator_DivBasisTimesScalar_impl_hpp__ 
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. 
 
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. 
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
EvaluatorStyle
An indication of how an Evaluator will behave. 
 
Teuchos::RCP< const PureBasis > getBasis() const 
 
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator. 
 
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...
 
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 ( ). 
 
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. 
 
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. 
 
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. 
 
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 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. 
 
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. 
 
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
 
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 ( ...