43#ifndef PANZER_DOF_GRADIENT_IMPL_HPP
44#define PANZER_DOF_GRADIENT_IMPL_HPP
50#include "Intrepid2_FunctionSpaceTools.hpp"
56 template <
typename ScalarT,
typename ArrayT>
57 struct evaluateGrad_withSens {
58 using scratch_view = Kokkos::View<ScalarT*,typename PHX::DevLayout<ScalarT>::type,
typename PHX::exec_space::scratch_memory_space,Kokkos::MemoryUnmanaged>;
59 using team_policy = Kokkos::TeamPolicy<PHX::exec_space>::member_type;
61 PHX::MDField<ScalarT> dof_grad_;
62 PHX::MDField<const ScalarT,Cell,Point> dof_value_;
63 const ArrayT grad_basis_;
64 const int num_fields_;
65 const int num_points_;
68 const bool use_shared_memory_;
70 evaluateGrad_withSens(PHX::MDField<ScalarT> dof_grad,
71 PHX::MDField<const ScalarT,Cell,Point> dof_value,
72 const ArrayT grad_basis,
73 const bool use_shared_memory) :
75 dof_value_(dof_value),
76 grad_basis_(grad_basis),
77 num_fields_(grad_basis.extent(1)),
78 num_points_(grad_basis.extent(2)),
79 space_dim_(grad_basis.extent(3)),
80 fad_size_(Kokkos::dimension_scalar(dof_grad.get_view())),
81 use_shared_memory_(use_shared_memory) {}
83 KOKKOS_INLINE_FUNCTION
84 void operator() (
const team_policy& team)
const
86 const int cell = team.league_rank();
88 if (not use_shared_memory_) {
89 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points_), [&] (
const int& pt) {
90 for (
int d=0; d<space_dim_; ++d) {
93 dof_grad_(cell,pt,d) = dof_value_(cell, 0) * grad_basis_(cell, 0, pt, d);
94 for (
int bf=1; bf<num_fields_; ++bf)
95 dof_grad_(cell,pt,d) += dof_value_(cell, bf) * grad_basis_(cell, bf, pt, d);
99 scratch_view dof_values;
100 scratch_view point_values;
101 if (Sacado::IsADType<ScalarT>::value) {
102 dof_values = scratch_view(team.team_shmem(),num_fields_,fad_size_);
103 point_values = scratch_view(team.team_shmem(),num_points_,fad_size_);
106 dof_values = scratch_view(team.team_shmem(),num_fields_);
107 point_values = scratch_view(team.team_shmem(),num_points_);
110 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_fields_), [&] (
const int& dof) {
111 dof_values(dof) = dof_value_(cell,dof);
115 for (
int dim=0; dim < space_dim_; ++dim) {
116 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points_), [&] (
const int& pt) {
117 point_values(pt) = 0.0;
120 for (
int dof=0; dof<num_fields_; ++dof) {
121 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points_), [&] (
const int& pt) {
122 point_values(pt) += dof_values(dof) * grad_basis_(cell,dof,pt,dim);
126 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points_), [&] (
const int& pt) {
127 dof_grad_(cell,pt,dim) = point_values(pt);
133 size_t team_shmem_size(
int )
const
135 if (not use_shared_memory_)
139 if (Sacado::IsADType<ScalarT>::value)
140 bytes = scratch_view::shmem_size(num_fields_,fad_size_) + scratch_view::shmem_size(num_points_,fad_size_);
142 bytes = scratch_view::shmem_size(num_fields_) + scratch_view::shmem_size(num_points_);
150template<
typename EvalT,
typename Traits>
153 const Teuchos::ParameterList& p) :
161 Teuchos::RCP<const PureBasis> basis
162 = p.get< Teuchos::RCP<BasisIRLayout> >(
"Basis")->getBasis();
165 TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsGrad(),std::logic_error,
166 "DOFGradient: Basis of type \"" << basis->name() <<
"\" does not support GRAD");
171 std::string n =
"DOFGradient: " +
dof_gradient.fieldTag().name() +
" ("+PHX::print<EvalT>()+
")";
176template<
typename EvalT,
typename TRAITS>
179 const PHX::FieldTag & output,
189 TEUCHOS_TEST_FOR_EXCEPTION(
bd_.getType()==
"HGrad",std::logic_error,
190 "DOFGradient: Basis of type \"" <<
bd_.getType() <<
"\" does not support GRAD");
195 std::string n =
"DOFGradient: " +
dof_gradient.fieldTag().name() +
" ("+PHX::print<EvalT>()+
")";
200template<
typename EvalT,
typename Traits>
212template<
typename EvalT,
typename Traits>
229 Kokkos::parallel_for(
"panzer::DOFGradient::evaluateFields", policy, eval);
ConstArray_CellBasisIPDim getGradBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the gradient of the basis evaluated at mesh points.
PHX::MDField< const Scalar, Cell, BASIS, IP, Dim > ConstArray_CellBasisIPDim
Array_CellBasisIPDim grad_basis
PHX::MDField< ScalarT > dof_gradient
panzer::IntegrationDescriptor id_
panzer::BasisDescriptor bd_
DOFGradient(const Teuchos::ParameterList &p)
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
PHX::MDField< const ScalarT, Cell, Point > dof_value
void evaluateFields(typename TRAITS::EvalData d)
WorksetDetailsAccessor wda
static HP & inst()
Private ctor.
bool useSharedMemory() const
Kokkos::TeamPolicy< TeamPolicyProperties... > teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
int num_cells
DEPRECATED - use: numCells()
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.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_
const panzer::Workset & EvalData