Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_DOFGradient_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef PANZER_DOF_GRADIENT_IMPL_HPP
12 #define PANZER_DOF_GRADIENT_IMPL_HPP
13 
15 #include "Panzer_BasisIRLayout.hpp"
18 #include "Intrepid2_FunctionSpaceTools.hpp"
19 
20 namespace panzer {
21 
22 namespace {
23 
24  template <typename ScalarT,typename ArrayT>
25  struct evaluateGrad_withSens {
26  using scratch_view = Kokkos::View<ScalarT*,typename PHX::DevLayout<ScalarT>::type,typename PHX::exec_space::scratch_memory_space,Kokkos::MemoryUnmanaged>;
27  using team_policy = Kokkos::TeamPolicy<PHX::exec_space>::member_type;
28 
31  const ArrayT grad_basis_;
32  const int num_fields_;
33  const int num_points_;
34  const int space_dim_;
35  const int fad_size_;
36  const bool use_shared_memory_;
37 
38  evaluateGrad_withSens(PHX::MDField<ScalarT> dof_grad,
40  const ArrayT grad_basis,
41  const bool use_shared_memory) :
42  dof_grad_(dof_grad),
43  dof_value_(dof_value),
44  grad_basis_(grad_basis),
45  num_fields_(grad_basis.extent(1)),
46  num_points_(grad_basis.extent(2)),
47  space_dim_(grad_basis.extent(3)),
48  fad_size_(Kokkos::dimension_scalar(dof_grad.get_view())),
49  use_shared_memory_(use_shared_memory) {}
50 
51  KOKKOS_INLINE_FUNCTION
52  void operator() (const team_policy& team) const
53  {
54  const int cell = team.league_rank();
55 
56  if (not use_shared_memory_) {
57  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points_), [&] (const int& pt) {
58  for (int d=0; d<space_dim_; ++d) {
59  // first initialize to the right thing (prevents over writing with 0)
60  // then loop over one less basis function
61  dof_grad_(cell,pt,d) = dof_value_(cell, 0) * grad_basis_(cell, 0, pt, d);
62  for (int bf=1; bf<num_fields_; ++bf)
63  dof_grad_(cell,pt,d) += dof_value_(cell, bf) * grad_basis_(cell, bf, pt, d);
64  }
65  });
66  } else {
67  scratch_view dof_values;
68  scratch_view point_values;
69  if (Sacado::IsADType<ScalarT>::value) {
70  dof_values = scratch_view(team.team_shmem(),num_fields_,fad_size_);
71  point_values = scratch_view(team.team_shmem(),num_points_,fad_size_);
72  }
73  else {
74  dof_values = scratch_view(team.team_shmem(),num_fields_);
75  point_values = scratch_view(team.team_shmem(),num_points_);
76  }
77 
78  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_fields_), [&] (const int& dof) {
79  dof_values(dof) = dof_value_(cell,dof);
80  });
81  team.team_barrier();
82 
83  for (int dim=0; dim < space_dim_; ++dim) {
84  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points_), [&] (const int& pt) {
85  point_values(pt) = 0.0;
86  });
87  // Perform contraction
88  for (int dof=0; dof<num_fields_; ++dof) {
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);
91  });
92  }
93  // Copy to main memory
94  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_points_), [&] (const int& pt) {
95  dof_grad_(cell,pt,dim) = point_values(pt);
96  });
97  } // loop over dim
98  }
99  }
100 
101  size_t team_shmem_size(int /* team_size */ ) const
102  {
103  if (not use_shared_memory_)
104  return 0;
105 
106  size_t bytes;
107  if (Sacado::IsADType<ScalarT>::value)
108  bytes = scratch_view::shmem_size(num_fields_,fad_size_) + scratch_view::shmem_size(num_points_,fad_size_);
109  else
110  bytes = scratch_view::shmem_size(num_fields_) + scratch_view::shmem_size(num_points_);
111  return bytes;
112  }
113  };
114 
115 } // anonymous namespace
116 
117 //**********************************************************************
118 template<typename EvalT, typename Traits>
121  const Teuchos::ParameterList& p) :
122  use_descriptors_(false),
123  dof_value( p.get<std::string>("Name"),
124  p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
125  dof_gradient( p.get<std::string>("Gradient Name"),
126  p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_vector ),
127  basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
128 {
130  = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
131 
132  // Verify that this basis supports the gradient operation
133  TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsGrad(),std::logic_error,
134  "DOFGradient: Basis of type \"" << basis->name() << "\" does not support GRAD");
135 
136  this->addEvaluatedField(dof_gradient);
137  this->addDependentField(dof_value);
138 
139  std::string n = "DOFGradient: " + dof_gradient.fieldTag().name() + " ("+PHX::print<EvalT>()+")";
140  this->setName(n);
141 }
142 
143 //**********************************************************************
144 template<typename EvalT, typename TRAITS>
147  const PHX::FieldTag & output,
148  const panzer::BasisDescriptor & bd,
150  : use_descriptors_(true)
151  , bd_(bd)
152  , id_(id)
153  , dof_value(input)
154  , dof_gradient(output)
155 {
156  // Verify that this basis supports the gradient operation
157  TEUCHOS_TEST_FOR_EXCEPTION(bd_.getType()=="HGrad",std::logic_error,
158  "DOFGradient: Basis of type \"" << bd_.getType() << "\" does not support GRAD");
159 
160  this->addEvaluatedField(dof_gradient);
161  this->addDependentField(dof_value);
162 
163  std::string n = "DOFGradient: " + dof_gradient.fieldTag().name() + " ("+PHX::print<EvalT>()+")";
164  this->setName(n);
165 }
166 
167 //**********************************************************************
168 template<typename EvalT, typename Traits>
169 void
172  typename Traits::SetupData sd,
173  PHX::FieldManager<Traits>& /* fm */)
174 {
175  if(not use_descriptors_)
176  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
177 }
178 
179 //**********************************************************************
180 template<typename EvalT, typename Traits>
181 void
184 {
185  if (workset.num_cells == 0)
186  return;
187 
188  const panzer::BasisValues2<double> & basisValues = use_descriptors_ ? this->wda(workset).getBasisValues(bd_,id_)
189  : *this->wda(workset).bases[basis_index];
190 
192  Array grad_basis = use_descriptors_ ? basisValues.getGradBasisValues(false) : Array(basisValues.grad_basis);
193 
194  bool use_shared_memory = panzer::HP::inst().useSharedMemory<ScalarT>();
195  evaluateGrad_withSens<ScalarT, Array> eval(dof_gradient,dof_value,grad_basis,use_shared_memory);
196  auto policy = panzer::HP::inst().teamPolicy<ScalarT,PHX::Device>(workset.num_cells);
197  Kokkos::parallel_for("panzer::DOFGradient::evaluateFields", policy, eval);
198 }
199 
200 //**********************************************************************
201 
202 }
203 
204 #endif
const ArrayT grad_basis_
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.
const int fad_size_
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.
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)
const int num_fields_
const int num_points_
PHX::MDField< const ScalarT, Cell, Point > dof_value_
const int space_dim_
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_