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;
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) :
80 fad_size_(Kokkos::dimension_scalar(dof_grad.get_view())),
83 KOKKOS_INLINE_FUNCTION
84 void operator() (
const team_policy& team)
const
86 const int cell = team.league_rank();
89 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_points_), [&] (
const int& pt) {
99 scratch_view dof_values;
100 scratch_view point_values;
101 if (Sacado::IsADType<ScalarT>::value) {
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) {
116 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_points_), [&] (
const int& pt) {
117 point_values(pt) = 0.0;
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)
142 bytes = scratch_view::shmem_size(num_fields_) + scratch_view::shmem_size(
num_points_);
150 template<
typename EvalT,
typename Traits>
154 use_descriptors_(false),
157 dof_gradient( p.get<std::string>(
"Gradient Name"),
166 "DOFGradient: Basis of type \"" << basis->name() <<
"\" does not support GRAD");
171 std::string n =
"DOFGradient: " +
dof_gradient.fieldTag().name() +
" ("+PHX::typeAsString<EvalT>()+
")";
176 template<
typename EvalT,
typename TRAITS>
179 const PHX::FieldTag & output,
182 : use_descriptors_(true)
186 , dof_gradient(output)
190 "DOFGradient: Basis of type \"" <<
bd_.
getType() <<
"\" does not support GRAD");
195 std::string n =
"DOFGradient: " +
dof_gradient.fieldTag().name() +
" ("+PHX::typeAsString<EvalT>()+
")";
200 template<
typename EvalT,
typename Traits>
207 if(not use_descriptors_)
212 template<
typename EvalT,
typename Traits>
221 : *this->wda(workset).bases[basis_index];
223 typedef decltype(basisValues.
grad_basis) ArrayT;
225 evaluateGrad_withSens<ScalarT, ArrayT> eval(dof_gradient,
dof_value,basisValues.
grad_basis,use_shared_memory);
227 Kokkos::parallel_for(policy, eval,
"panzer::DOFGradient::evaluateFields");
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
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
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_