11#ifndef PANZER_DOF_IMPL_HPP
12#define PANZER_DOF_IMPL_HPP
22#include "Intrepid2_FunctionSpaceTools.hpp"
35template<
typename EvalT,
typename TRAITS>
39 const std::string fieldName = p.get<std::string>(
"Name");
40 basis = p.get< Teuchos::RCP<const PureBasis> >(
"Basis");
41 Teuchos::RCP<const PointRule> pointRule = p.get< Teuchos::RCP<const PointRule> >(
"Point Rule");
42 is_vector_basis = basis->isVectorBasis();
44 std::string evalName = fieldName+
"_"+pointRule->getName();
45 if(p.isType<
bool>(
"Use DOF Name")) {
46 if(p.get<
bool>(
"Use DOF Name"))
50 dof_basis = PHX::MDField<const ScalarT,Cell,Point>(fieldName, basis->functional);
52 this->addDependentField(dof_basis);
55 Teuchos::RCP<BasisIRLayout> layout = Teuchos::rcp(
new BasisIRLayout(basis,*pointRule));
56 basisValues = Teuchos::rcp(
new BasisValues2<double>(basis->name()+
"_"+pointRule->getName()+
"_"));
57 basisValues->setupArrays(layout,
false);
61 if(basis->isScalarBasis()) {
62 dof_ip_scalar = PHX::MDField<ScalarT,Cell,Point>(
64 pointRule->dl_scalar);
65 this->addEvaluatedField(dof_ip_scalar);
66 this->addNonConstDependentField(basisValues->basis_ref_scalar);
67 this->addNonConstDependentField(basisValues->basis_scalar);
69 else if(basis->isVectorBasis()) {
70 dof_ip_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(
72 pointRule->dl_vector);
73 this->addEvaluatedField(dof_ip_vector);
74 this->addNonConstDependentField(basisValues->basis_ref_vector);
75 this->addNonConstDependentField(basisValues->basis_vector);
78 { TEUCHOS_ASSERT(
false); }
80 std::string n =
"DOF_PointValues: " + dof_basis.fieldTag().name();
85template<
typename EvalT,
typename TRAITS>
90 if(!is_vector_basis) {
91 this->utils.setFieldData(basisValues->basis_ref_scalar,fm);
92 this->utils.setFieldData(basisValues->basis_scalar,fm);
95 this->utils.setFieldData(basisValues->basis_ref_vector,fm);
96 this->utils.setFieldData(basisValues->basis_vector,fm);
101template<
typename EvalT,
typename TRAITS>
107 if(is_vector_basis) {
108 int spaceDim = basisValues->basis_vector.extent(3);
111 Kokkos::parallel_for(Kokkos::TeamPolicy<PHX::Device>(workset.num_cells,Kokkos::AUTO(),vector_size),functor);
115 Kokkos::parallel_for(Kokkos::TeamPolicy<PHX::Device>(workset.num_cells,Kokkos::AUTO(),vector_size),functor);
120 Kokkos::parallel_for(workset.num_cells,functor);
131template<
typename TRAITS>
135 const std::string fieldName = p.get<std::string>(
"Name");
136 basis = p.get< Teuchos::RCP<const PureBasis> >(
"Basis");
137 Teuchos::RCP<const PointRule> pointRule = p.get< Teuchos::RCP<const PointRule> >(
"Point Rule");
138 is_vector_basis = basis->isVectorBasis();
140 if(p.isType<Teuchos::RCP<
const std::vector<int> > >(
"Jacobian Offsets Vector")) {
141 const std::vector<int> &
offsets = *p.get<Teuchos::RCP<const std::vector<int> > >(
"Jacobian Offsets Vector");
144 offsets_array = PHX::View<int*>(
"offsets",
offsets.size());
145 for(std::size_t i=0;i<
offsets.size();i++)
148 accelerate_jacobian =
true;
151 accelerate_jacobian =
false;
153 std::string evalName = fieldName+
"_"+pointRule->getName();
154 if(p.isType<
bool>(
"Use DOF Name")) {
155 if(p.get<
bool>(
"Use DOF Name"))
156 evalName = fieldName;
159 dof_basis = PHX::MDField<const ScalarT,Cell,Point>(fieldName, basis->functional);
161 this->addDependentField(dof_basis);
164 Teuchos::RCP<BasisIRLayout> layout = Teuchos::rcp(
new BasisIRLayout(basis,*pointRule));
165 basisValues = Teuchos::rcp(
new BasisValues2<double>(basis->name()+
"_"+pointRule->getName()+
"_"));
166 basisValues->setupArrays(layout,
false);
170 if(basis->isScalarBasis()) {
171 dof_ip_scalar = PHX::MDField<ScalarT,Cell,Point>(
173 pointRule->dl_scalar);
174 this->addEvaluatedField(dof_ip_scalar);
175 this->addNonConstDependentField(basisValues->basis_ref_scalar);
176 this->addNonConstDependentField(basisValues->basis_scalar);
178 else if(basis->isVectorBasis()) {
179 dof_ip_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(
181 pointRule->dl_vector);
182 this->addEvaluatedField(dof_ip_vector);
183 this->addNonConstDependentField(basisValues->basis_ref_vector);
184 this->addNonConstDependentField(basisValues->basis_vector);
187 { TEUCHOS_ASSERT(
false); }
189 std::string n =
"DOF_PointValues: " + dof_basis.fieldTag().name() +
" Jacobian";
194template<
typename TRAITS>
199 if(!is_vector_basis) {
200 this->utils.setFieldData(basisValues->basis_ref_scalar,fm);
201 this->utils.setFieldData(basisValues->basis_scalar,fm);
204 this->utils.setFieldData(basisValues->basis_ref_vector,fm);
205 this->utils.setFieldData(basisValues->basis_vector,fm);
210template<
typename TRAITS>
216 if(is_vector_basis) {
217 if(accelerate_jacobian) {
218 int spaceDim = basisValues->basis_vector.extent(3);
220 dof_functors::EvaluateDOFFastSens_Vector<ScalarT,typename BasisValues2<double>::Array_CellBasisIPDim,3> functor(dof_basis,dof_ip_vector,offsets_array,basisValues->basis_vector);
221 Kokkos::parallel_for(workset.num_cells,functor);
224 dof_functors::EvaluateDOFFastSens_Vector<ScalarT,typename BasisValues2<double>::Array_CellBasisIPDim,2> functor(dof_basis,dof_ip_vector,offsets_array,basisValues->basis_vector);
225 Kokkos::parallel_for(workset.num_cells,functor);
229 int spaceDim = basisValues->basis_vector.extent(3);
231 dof_functors::EvaluateDOFWithSens_Vector<ScalarT,typename BasisValues2<double>::Array_CellBasisIPDim,3> functor(dof_basis.get_static_view(),dof_ip_vector.get_static_view(),basisValues->basis_vector);
232 Kokkos::parallel_for(Kokkos::TeamPolicy<PHX::Device>(workset.num_cells,Kokkos::AUTO(),vector_size),functor);
235 dof_functors::EvaluateDOFWithSens_Vector<ScalarT,typename BasisValues2<double>::Array_CellBasisIPDim,2> functor(dof_basis.get_static_view(),dof_ip_vector.get_static_view(),basisValues->basis_vector);
236 Kokkos::parallel_for(Kokkos::TeamPolicy<PHX::Device>(workset.num_cells,Kokkos::AUTO(),vector_size),functor);
241 if(accelerate_jacobian) {
242 dof_functors::EvaluateDOFFastSens_Scalar<ScalarT,typename BasisValues2<double>::Array_CellBasisIP> functor(dof_basis,dof_ip_scalar,offsets_array,basisValues->basis_scalar);
243 Kokkos::parallel_for(workset.num_cells,functor);
246 dof_functors::EvaluateDOFWithSens_Scalar<ScalarT,typename BasisValues2<double>::Array_CellBasisIP> functor(dof_basis,dof_ip_scalar,basisValues->basis_scalar);
247 Kokkos::parallel_for(workset.num_cells,functor);
PHX::View< const int * > offsets
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
void evaluateFields(typename TRAITS::EvalData d)
DOF_PointValues(const Teuchos::ParameterList &p)
int vectorSize() const
Returns the vector size. Specialized for AD scalar types.
static HP & inst()
Private ctor.