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())
 
   87     using PHX::DataLayout;
 
   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_);
 
  129       View<View<const ScalarT**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>(
"GradBasisDotVector::KokkosFieldMultipliers",
 
  131     for (
const auto& name : fmNames)
 
  138     string n(
"Integrator_GradBasisDotVector (");
 
  143     n += 
"):  " + 
field_.fieldTag().name();
 
  152   template<
typename EvalT, 
typename Traits>
 
  159       p.get<std::string>(
"Residual Name"),
 
  160       p.get<std::string>(
"Flux 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>(),
 
  168       p.isType<Teuchos::
RCP<PHX::DataLayout>>(
"Vector Data Layout") ?
 
  169         p.get<Teuchos::
RCP<PHX::DataLayout>>(
"Vector Data Layout") :
 
  185   template<
typename EvalT, 
typename Traits>
 
  197     for (
size_t i(0); i < fieldMults_.size(); ++i)
 
  198       kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
 
  210   template<
typename EvalT, 
typename Traits>
 
  211   template<
int NUM_FIELD_MULT>
 
  212   KOKKOS_INLINE_FUNCTION
 
  217     const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
 const 
  220     const int cell = team.league_rank();
 
  223     const int numQP(vector_.extent(1)), numDim(vector_.extent(2)),
 
  224               numBases(basis_.extent(1));
 
  226       Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  227         field_(cell, basis) = 0.0;
 
  233     if (NUM_FIELD_MULT == 0)
 
  238       for (
int qp(0); qp < numQP; ++qp)
 
  240         for (
int dim(0); dim < numDim; ++dim)
 
  242           tmp = multiplier_ * vector_(cell, qp, dim);
 
  243     Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  244             field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
 
  249     else if (NUM_FIELD_MULT == 1)
 
  254       for (
int qp(0); qp < numQP; ++qp)
 
  256         for (
int dim(0); dim < numDim; ++dim)
 
  258           tmp = multiplier_ * vector_(cell, qp, dim) *
 
  259             kokkosFieldMults_(0)(cell, qp);
 
  260     Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  261             field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
 
  273       const int numFieldMults(kokkosFieldMults_.extent(0));
 
  274       for (
int qp(0); qp < numQP; ++qp)
 
  277         for (
int fm(0); fm < numFieldMults; ++fm)
 
  278           fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
 
  279         for (
int dim(0); dim < numDim; ++dim)
 
  281           tmp = multiplier_ * vector_(cell, qp, dim) * fieldMultsTotal;
 
  282     Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  283             field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
 
  295   template<
typename EvalT, 
typename Traits>
 
  296   template<
int NUM_FIELD_MULT>
 
  297   KOKKOS_INLINE_FUNCTION
 
  302     const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
 const 
  305     const int cell = team.league_rank();
 
  306     const int numQP = vector_.extent(1);
 
  307     const int numDim = vector_.extent(2);
 
  308     const int numBases = basis_.extent(1);
 
  309     const int fadSize = Kokkos::dimension_scalar(field_.get_view());
 
  313     if (Sacado::IsADType<ScalarT>::value) {
 
  315       tmp_field = 
scratch_view(team.team_shmem(),numBases,fadSize);
 
  323     Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  324       tmp_field(basis) = 0.0;
 
  329     if (NUM_FIELD_MULT == 0)
 
  334       for (
int qp(0); qp < numQP; ++qp)
 
  336         for (
int dim(0); dim < numDim; ++dim)
 
  339           tmp(0) = multiplier_ * vector_(cell, qp, dim);
 
  340     Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  341       tmp_field(basis) += basis_(cell, basis, qp, dim) * tmp(0);
 
  346     else if (NUM_FIELD_MULT == 1)
 
  351       for (
int qp(0); qp < numQP; ++qp)
 
  353         for (
int dim(0); dim < numDim; ++dim)
 
  356           tmp(0) = multiplier_ * vector_(cell, qp, dim) *
 
  357             kokkosFieldMults_(0)(cell, qp);
 
  358     Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  359         tmp_field(basis) += basis_(cell, basis, qp, dim) * tmp(0);
 
  371       const int numFieldMults(kokkosFieldMults_.extent(0));
 
  372       for (
int qp(0); qp < numQP; ++qp)
 
  375         for (
int fm(0); fm < numFieldMults; ++fm)
 
  376           fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
 
  377         for (
int dim(0); dim < numDim; ++dim)
 
  380           tmp(0) = multiplier_ * vector_(cell, qp, dim) * fieldMultsTotal;
 
  381     Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (
const int& basis) {
 
  382       tmp_field(basis) += basis_(cell, basis, qp, dim) * tmp(0);
 
  390       Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int& basis) {
 
  391   field_(cell,basis) = tmp_field(basis);
 
  395       Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int& basis) {
 
  396   field_(cell,basis) += tmp_field(basis);
 
  407   template<
typename EvalT, 
typename Traits>
 
  413     using Kokkos::parallel_for;
 
  414     using Kokkos::TeamPolicy;
 
  417     basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
 
  420     if (use_shared_memory) {
 
  422       if (Sacado::IsADType<ScalarT>::value) {
 
  423   const int fadSize = Kokkos::dimension_scalar(field_.get_view());
 
  424   bytes = scratch_view::shmem_size(1,fadSize) + scratch_view::shmem_size(basis_.extent(1),fadSize);
 
  427   bytes = scratch_view::shmem_size(1) + scratch_view::shmem_size(basis_.extent(1));
 
  432       if (fieldMults_.size() == 0) {
 
  434         parallel_for(policy, *
this, this->getName());
 
  435       } 
else if (fieldMults_.size() == 1) {
 
  437         parallel_for(policy, *
this, this->getName());
 
  440         parallel_for(policy, *
this, this->getName());
 
  447       if (fieldMults_.size() == 0) {
 
  449         parallel_for(policy, *
this, this->getName());
 
  450       } 
else if (fieldMults_.size() == 1) {
 
  452         parallel_for(policy, *
this, this->getName());
 
  455         parallel_for(policy, *
this, this->getName());
 
  465   template<
typename EvalT, 
typename TRAITS>
 
  472     using PHX::DataLayout;
 
  481     p->set<
string>(
"Residual Name", 
"?");
 
  482     p->set<
string>(
"Flux Name", 
"?");
 
  484     p->set(
"Basis", basis);
 
  487     p->set<
double>(
"Multiplier", 1.0);
 
  489     p->set(
"Field Multipliers", fms);
 
  491     p->set(
"Vector Data Layout", vecDL);
 
  497 #endif // __Panzer_Integrator_GradBasisDotVector_impl_hpp__ 
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. 
 
Kokkos::View< ScalarT *,typename PHX::DevLayout< ScalarT >::type, typename PHX::exec_space::scratch_memory_space, Kokkos::MemoryUnmanaged > scratch_view
Type for shared memory. 
 
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. 
 
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> 
 
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...
 
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_