11 #ifndef PANZER_DOF_DIV_IMPL_HPP
12 #define PANZER_DOF_DIV_IMPL_HPP
22 template<
typename ScalarT,
typename ArrayT>
37 const ArrayT & in_div_basis,
38 const bool in_use_shared_memory =
false)
44 fad_size(static_cast<int>(Kokkos::dimension_scalar(in_dof_div.get_view()))),
48 KOKKOS_INLINE_FUNCTION
49 void operator()(
const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
51 const int cell = team.league_rank();
57 if (Sacado::IsADType<ScalarT>::value) {
66 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_fields), [&] (
const int& dof) {
70 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_points), [&] (
const int& pt) {
71 point_values(pt) = 0.0;
78 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_points), [&] (
const int& pt) {
79 point_values(pt) += dof_values(dof) *
div_basis(cell,dof,pt);
84 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_points), [&] (
const int& pt) {
85 dof_div(cell,pt) = point_values(pt);
89 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,
num_points), [&] (
const int& pt) {
105 if (Sacado::IsADType<ScalarT>::value)
119 template<
typename EvalT,
typename TRAITS>
122 use_descriptors_(false),
132 "DOFDiv: Basis of type \"" << basis->name() <<
"\" does not support DIV");
134 "DOFDiv: Basis of type \"" << basis->name() <<
"\" in DOF Div should require orientations. So we are throwing.");
141 this->addEvaluatedField(dof_div);
144 std::string n =
"DOFDiv: " + dof_div.fieldTag().name() +
" ("+PHX::print<EvalT>()+
")";
149 template<
typename EvalT,
typename TRAITS>
155 : use_descriptors_(true)
166 this->addEvaluatedField(dof_div);
169 std::string n =
"DOFDiv: " + dof_div.fieldTag().name() +
" ("+PHX::print<EvalT>()+
")";
174 template<
typename EvalT,
typename TRAITS>
179 if(not use_descriptors_)
184 template<
typename EvalT,
typename TRAITS>
188 if (workset.num_cells == 0)
192 : *this->wda(workset).bases[basis_index];
200 Kokkos::parallel_for(this->getName(),policy,f);
210 template<
typename TRAITS>
213 use_descriptors_(false),
225 accelerate_jacobian =
true;
228 accelerate_jacobian =
false;
229 accelerate_jacobian =
false;
233 "DOFDiv: Basis of type \"" << basis->name() <<
"\" does not support DIV");
235 "DOFDiv: Basis of type \"" << basis->name() <<
"\" in DOF Div should require orientations. So we are throwing.");
242 this->addEvaluatedField(dof_div);
245 std::string n =
"DOFDiv: " + dof_div.fieldTag().name() +
" ("+PHX::print<panzer::Traits::Jacobian>()+
")";
250 template<
typename TRAITS>
256 : use_descriptors_(true)
266 accelerate_jacobian =
false;
269 this->addEvaluatedField(dof_div);
272 std::string n =
"DOFDiv: " + dof_div.fieldTag().name() +
" ("+PHX::print<panzer::Traits::Jacobian>()+
")";
277 template<
typename TRAITS>
283 this->utils.setFieldData(
dof_div,fm);
289 template<
typename TRAITS>
293 if (workset.num_cells == 0)
302 if(!accelerate_jacobian) {
306 Kokkos::parallel_for(this->
getName(),policy,f);
panzer::BasisDescriptor bd_
panzer::IntegrationDescriptor id_
DOFDiv(const Teuchos::ParameterList &p)
T & get(const std::string &name, T def_value)
Kokkos::TeamPolicy< TeamPolicyProperties...> teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
Array_CellBasisIP div_basis
PHX::View< const int * > offsets
panzer::Traits::Jacobian::ScalarT ScalarT
PHX::MDField< ScalarT, Cell, IP > dof_div
const std::string & getType() const
Get type of basis.
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
virtual const std::string & getName() const =0
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.
ConstArray_CellBasisIP getDivVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the divergence of a vector basis evaluated at mesh points.
KOKKOS_INLINE_FUNCTION void operator()(const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
void evaluateFields(typename TRAITS::EvalData d)
PHX::MDField< const ScalarT, Cell, Point > dof_value
bool isType(const std::string &name) const
WorksetDetailsAccessor wda
static HP & inst()
Private ctor.
#define TEUCHOS_ASSERT(assertion_test)
PHX::MDField< const ScalarT, Cell, Point > dof_value
size_t team_shmem_size(int) const
EvaluateDOFDiv_withSens(PHX::MDField< ScalarT, Cell, IP > &in_dof_div, PHX::MDField< const ScalarT, Cell, Point > &in_dof_value, const ArrayT &in_div_basis, const bool in_use_shared_memory=false)
const bool use_shared_memory
PHX::MDField< ScalarT, Cell, IP > dof_div
Kokkos::View< ScalarT *,typename PHX::DevLayout< ScalarT >::type, typename PHX::exec_space::scratch_memory_space, Kokkos::MemoryUnmanaged > scratch_view