11 #ifndef PANZER_DOF_GRADIENT_IMPL_HPP
12 #define PANZER_DOF_GRADIENT_IMPL_HPP
18 #include "Intrepid2_FunctionSpaceTools.hpp"
24 template <
typename ScalarT,
typename ArrayT>
25 struct evaluateGrad_withSens {
27 using team_policy = Kokkos::TeamPolicy<PHX::exec_space>::member_type;
40 const ArrayT grad_basis,
41 const bool use_shared_memory) :
48 fad_size_(Kokkos::dimension_scalar(dof_grad.get_view())),
51 KOKKOS_INLINE_FUNCTION
52 void operator() (
const team_policy& team)
const
54 const int cell = team.league_rank();
57 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_points_), [&] (
const int& pt) {
67 scratch_view dof_values;
68 scratch_view point_values;
69 if (Sacado::IsADType<ScalarT>::value) {
74 dof_values = scratch_view(team.team_shmem(),
num_fields_);
75 point_values = scratch_view(team.team_shmem(),
num_points_);
78 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_fields_), [&] (
const int& dof) {
84 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_points_), [&] (
const int& pt) {
85 point_values(pt) = 0.0;
89 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_points_), [&] (
const int& pt) {
90 point_values(pt) += dof_values(dof) *
grad_basis_(cell,dof,pt,dim);
94 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_points_), [&] (
const int& pt) {
95 dof_grad_(cell,pt,dim) = point_values(pt);
101 size_t team_shmem_size(
int )
const
103 if (not use_shared_memory_)
107 if (Sacado::IsADType<ScalarT>::value)
110 bytes = scratch_view::shmem_size(num_fields_) + scratch_view::shmem_size(
num_points_);
118 template<
typename EvalT,
typename Traits>
122 use_descriptors_(false),
125 dof_gradient( p.get<std::string>(
"Gradient Name"),
134 "DOFGradient: Basis of type \"" << basis->name() <<
"\" does not support GRAD");
139 std::string n =
"DOFGradient: " +
dof_gradient.fieldTag().name() +
" ("+PHX::print<EvalT>()+
")";
144 template<
typename EvalT,
typename TRAITS>
150 : use_descriptors_(true)
154 , dof_gradient(output)
158 "DOFGradient: Basis of type \"" <<
bd_.
getType() <<
"\" does not support GRAD");
163 std::string n =
"DOFGradient: " +
dof_gradient.fieldTag().name() +
" ("+PHX::print<EvalT>()+
")";
168 template<
typename EvalT,
typename Traits>
175 if(not use_descriptors_)
180 template<
typename EvalT,
typename Traits>
189 : *this->wda(workset).bases[basis_index];
195 evaluateGrad_withSens<ScalarT, Array> eval(dof_gradient,
dof_value,grad_basis,use_shared_memory);
197 Kokkos::parallel_for(
"panzer::DOFGradient::evaluateFields", policy, eval);
int num_cells
DEPRECATED - use: numCells()
const bool use_shared_memory_
T & get(const std::string &name, T def_value)
Kokkos::TeamPolicy< TeamPolicyProperties...> teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
DOFGradient(const Teuchos::ParameterList &p)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
PHX::MDField< ScalarT > dof_gradient
PHX::MDField< ScalarT > dof_grad_
panzer::BasisDescriptor bd_
const std::string & getType() const
Get type of basis.
PHX::MDField< const ScalarT, Cell, Point > dof_value
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
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)
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
PHX::MDField< const ScalarT, Cell, Point > dof_value_
Array_CellBasisIPDim grad_basis
static HP & inst()
Private ctor.
PHX::MDField< const ScalarT, Cell, Point > dof_value
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.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_