43 #ifndef PANZER_DOF_POINT_FIELD_DECL_HPP 44 #define PANZER_DOF_POINT_FIELD_DECL_HPP 49 #include "Intrepid2_FunctionSpaceTools.hpp" 51 #include "Phalanx_DataLayout_MDALayout.hpp" 56 template <
typename EvalT,
typename TRAITST>
59 const std::string & coordinateName,
60 const Teuchos::RCP<PHX::DataLayout> & coordLayout,
61 const Teuchos::RCP<PHX::DataLayout> & quadLayout,
62 const std::string & postfixFieldName)
66 int cellCount = fieldBasis.
functional->dimension(0);
67 int coeffCount = fieldBasis.
functional->dimension(1);
68 int pointCount = coordLayout->dimension(0);
69 int dimCount = coordLayout->dimension(1);
71 Teuchos::RCP<PHX::DataLayout> basisLayout = fieldBasis.
functional;
73 coordinates = PHX::MDField<ScalarT,Point,Dim>(coordinateName,coordLayout);
74 dof_coeff = PHX::MDField<ScalarT>(fieldName,basisLayout);
75 dof_field = PHX::MDField<ScalarT>(fieldName+postfixFieldName,quadLayout);
77 this->addDependentField(coordinates);
78 this->addDependentField(dof_coeff);
79 this->addEvaluatedField(dof_field);
82 basisRef = Kokkos::DynRankView<double,PHX::Device>(
"basisRef",coeffCount,pointCount);
83 basis = Kokkos::DynRankView<double,PHX::Device>(
"basis",cellCount,coeffCount,pointCount);
84 intrpCoords = Kokkos::DynRankView<double,PHX::Device>(
"intrpCoords",pointCount,dimCount);
86 std::string n =
"DOF_PointField: " + dof_field.fieldTag().name();
91 template <
typename EvalT,
typename TRAITST>
95 this->utils.setFieldData(coordinates,fm);
96 this->utils.setFieldData(dof_coeff,fm);
97 this->utils.setFieldData(dof_field,fm);
101 template <
typename EvalT,
typename TRAITST>
105 dof_field.deep_copy(
ScalarT(0.0));
108 for (
int i = 0; i < coordinates.extent_int(0); ++i)
109 for (
int j = 0; j < coordinates.extent_int(1); ++j)
110 intrpCoords(i,j) = Sacado::ScalarValue<ScalarT>::eval(coordinates(i,j));
112 if(workset.num_cells>0) {
114 intrepidBasis->getValues(basisRef, intrpCoords, Intrepid2::OPERATOR_VALUE);
117 Intrepid2::FunctionSpaceTools::
118 HGRADtransformVALUE<double>(
basis,
122 Intrepid2::FunctionSpaceTools::
123 evaluate<ScalarT>(dof_field,dof_coeff,
basis);
void initialize(const std::string &fieldName, const PureBasis &fieldBasis, const std::string &coordinateName, const Teuchos::RCP< PHX::DataLayout > &coordLayout, const Teuchos::RCP< PHX::DataLayout > &quadLayout, const std::string &postfixFieldName)
Convenience initialization routine, see constructor above.
void evaluateFields(typename TRAITST::EvalData workset)
Teuchos::RCP< Intrepid2::Basis< double, Kokkos::DynRankView< double, PHX::Device > > > getIntrepid2Basis() const
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
void postRegistrationSetup(typename TRAITST::SetupData d, PHX::FieldManager< TRAITST > &vm)
Description and data layouts associated with a particular basis.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis> or <Cell,Basis>