43 #ifndef PANZER_DOF_CURL_IMPL_HPP 
   44 #define PANZER_DOF_CURL_IMPL_HPP 
   49 #include "Intrepid2_FunctionSpaceTools.hpp" 
   50 #include "Phalanx_KokkosDeviceTypes.hpp" 
   57 template <
typename ScalarT,
typename Array,
int spaceDim>
 
   58 class EvaluateCurlWithSens_Vector {
 
   67   typedef typename PHX::Device execution_space;
 
   69   EvaluateCurlWithSens_Vector(PHX::MDField<const ScalarT,Cell,Point> in_dof_value,
 
   70                               PHX::MDField<ScalarT,Cell,Point,Dim> in_dof_curl,
 
   77   KOKKOS_INLINE_FUNCTION
 
   78   void operator()(
const unsigned int cell)
 const 
   81       for (
int d=0; d<spaceDim; d++) {
 
   92 template <
typename ScalarT,
typename ArrayT>
 
   93 void evaluateCurl_withSens_vector(
int numCells,
 
   94                            PHX::MDField<ScalarT,Cell,Point,Dim> & 
dof_curl, 
 
   95                            PHX::MDField<const ScalarT,Cell,Point> & 
dof_value,
 
  101     int numPoints = curl_basis.extent(2);
 
  102     int spaceDim  = curl_basis.extent(3);
 
  104     for (
int cell=0; cell<numCells; cell++) {
 
  106         for (
int d=0; d<spaceDim; d++) {
 
  119 template <
typename ScalarT,
typename Array>
 
  120 class EvaluateCurlWithSens_Scalar {
 
  121   PHX::MDField<const ScalarT,Cell,Point> 
dof_value;
 
  122   PHX::MDField<ScalarT,Cell,Point> 
dof_curl;
 
  129   typedef typename PHX::Device execution_space;
 
  131   EvaluateCurlWithSens_Scalar(PHX::MDField<const ScalarT,Cell,Point> in_dof_value,
 
  132                               PHX::MDField<ScalarT,Cell,Point> in_dof_curl,
 
  134     : dof_value(in_dof_value), dof_curl(in_dof_curl), curl_basis(in_curl_basis)
 
  136     numFields = curl_basis.extent(1);
 
  137     numPoints = curl_basis.extent(2);
 
  139   KOKKOS_INLINE_FUNCTION
 
  140   void operator()(
const unsigned int cell)
 const 
  152 template <
typename ScalarT,
typename ArrayT>
 
  153 void evaluateCurl_withSens_scalar(
int numCells,
 
  154                            PHX::MDField<ScalarT,Cell,Point> & dof_curl, 
 
  155                            PHX::MDField<const ScalarT,Cell,Point> & dof_value,
 
  156                            const ArrayT & curl_basis)
 
  160     int numFields = curl_basis.extent(1);
 
  161     int numPoints = curl_basis.extent(2);
 
  163     for (
int cell=0; cell<numCells; cell++) {
 
  176 template <
typename ScalarT,
typename Array,
int spaceDim>
 
  177 class EvaluateCurlFastSens_Vector {
 
  178   PHX::MDField<const ScalarT,Cell,Point> 
dof_value;
 
  179   PHX::MDField<ScalarT,Cell,Point,Dim> 
dof_curl;
 
  187   typedef typename PHX::Device execution_space;
 
  189   EvaluateCurlFastSens_Vector(PHX::MDField<const ScalarT,Cell,Point> in_dof_value,
 
  190                               PHX::MDField<ScalarT,Cell,Point,Dim> in_dof_curl,
 
  191                               Kokkos::View<const int*,PHX::Device> in_offsets,
 
  193     : dof_value(in_dof_value), dof_curl(in_dof_curl), 
offsets(in_offsets), curl_basis(in_curl_basis)
 
  195     numFields = curl_basis.extent(1);
 
  196     numPoints = curl_basis.extent(2);
 
  198   KOKKOS_INLINE_FUNCTION
 
  199   void operator()(
const unsigned int cell)
 const 
  202       for (
int d=0; d<spaceDim; d++) {
 
  215 template <
typename ScalarT,
typename ArrayT>
 
  216 void evaluateCurl_fastSens_vector(
int numCells,
 
  217                            PHX::MDField<ScalarT,Cell,Point,Dim> & dof_curl, 
 
  218                            PHX::MDField<const ScalarT,Cell,Point> & dof_value,
 
  219                            const std::vector<int> & 
offsets,
 
  220                            const ArrayT & curl_basis)
 
  223     int numFields = curl_basis.extent(1);
 
  224     int numPoints = curl_basis.extent(2);
 
  225     int spaceDim  = curl_basis.extent(3);
 
  227     for (
int cell=0; cell<numCells; cell++) {
 
  229         for (
int d=0; d<spaceDim; d++) {
 
  245 template <
typename ScalarT,
typename Array>
 
  246 class EvaluateCurlFastSens_Scalar {
 
  247   PHX::MDField<const ScalarT,Cell,Point> 
dof_value;
 
  248   PHX::MDField<ScalarT,Cell,Point> 
dof_curl;
 
  249   Kokkos::View<const int*,PHX::Device> 
offsets;
 
  256   typedef typename PHX::Device execution_space;
 
  258   EvaluateCurlFastSens_Scalar(PHX::MDField<const ScalarT,Cell,Point> in_dof_value,
 
  259                               PHX::MDField<ScalarT,Cell,Point> in_dof_curl,
 
  260                               Kokkos::View<const int*,PHX::Device> in_offsets,
 
  262     : dof_value(in_dof_value), dof_curl(in_dof_curl), offsets(in_offsets), curl_basis(in_curl_basis)
 
  264     numFields = curl_basis.extent(1);
 
  265     numPoints = curl_basis.extent(2);
 
  267   KOKKOS_INLINE_FUNCTION
 
  268   void operator()(
const unsigned int cell)
 const 
  282 template <
typename ScalarT,
typename ArrayT>
 
  283 void evaluateCurl_fastSens_scalar(
int numCells,
 
  284                            PHX::MDField<ScalarT,Cell,Point> & dof_curl, 
 
  285                            PHX::MDField<const ScalarT,Cell,Point> & dof_value,
 
  286                            const std::vector<int> & offsets,
 
  287                            const ArrayT & curl_basis)
 
  290     int numFields = curl_basis.extent(1);
 
  291     int numPoints = curl_basis.extent(2);
 
  293     for (
int cell=0; cell<numCells; cell++) {
 
  317 template<
typename EvalT, 
typename TRAITS>                   
 
  320   use_descriptors_(false),
 
  321   dof_value( p.get<std::string>(
"Name"), 
 
  330                              "DOFCurl: Basis of type \"" << basis->name() << 
"\" does not support CURL");
 
  332                              "DOFCurl: Basis of type \"" << basis->name() << 
"\" in DOF Curl should require orientations. So we are throwing.");
 
  336   if(basis_dimension==2) {
 
  341   else if(basis_dimension==3) {
 
  342      dof_curl_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(p.
get<std::string>(
"Curl Name"), 
 
  350   this->addDependentField(dof_value);
 
  357 template<
typename EvalT, 
typename TRAITS>                   
 
  360         const PHX::FieldTag & output,
 
  364   : use_descriptors_(true)
 
  374   if(basis_dimension==2) {
 
  378   else if(basis_dimension==3) {
 
  386   this->addDependentField(dof_value);
 
  393 template<
typename EvalT, 
typename TRAITS>                   
 
  398   this->utils.setFieldData(dof_value,fm);
 
  399   if(basis_dimension==3)
 
  400     this->utils.setFieldData(dof_curl_vector,fm);
 
  402     this->utils.setFieldData(dof_curl_scalar,fm);
 
  404   if(not use_descriptors_)
 
  409 template<
typename EvalT, 
typename TRAITS>                   
 
  414                                                                       : *this->wda(workset).bases[basis_index];
 
  416   if(basis_dimension==3) {
 
  417     EvaluateCurlWithSens_Vector<ScalarT,typename BasisValues2<double>::Array_CellBasisIPDim,3> functor(dof_value,dof_curl_vector,basisValues.
curl_basis_vector);
 
  418     Kokkos::parallel_for(workset.num_cells,functor);
 
  421     EvaluateCurlWithSens_Scalar<ScalarT,typename BasisValues2<double>::Array_CellBasisIP> functor(dof_value,dof_curl_scalar,basisValues.
curl_basis_scalar);
 
  422     Kokkos::parallel_for(workset.num_cells,functor);
 
  433 template<
typename TRAITS>                   
 
  436   use_descriptors_(false),
 
  437   dof_value( p.get<std::string>(
"Name"), 
 
  450     Kokkos::View<int*,PHX::Device> offsets_array_nc
 
  451         = Kokkos::View<int*,PHX::Device>(
"offsets",offsets.size());
 
  452     for(std::size_t i=0;i<offsets.size();i++)
 
  453       offsets_array_nc(i) = offsets[i];
 
  454     offsets_array = offsets_array_nc;
 
  456     accelerate_jacobian = 
true;  
 
  459     accelerate_jacobian = 
false; 
 
  463                              "DOFCurl: Basis of type \"" << basis->name() << 
"\" does not support CURL");
 
  465                              "DOFCurl: Basis of type \"" << basis->name() << 
"\" in DOF Curl should require orientations. So we are throwing.");
 
  469   if(basis_dimension==2) {
 
  474   else if(basis_dimension==3) {
 
  475      dof_curl_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(p.
get<std::string>(
"Curl Name"), 
 
  483   this->addDependentField(dof_value);
 
  490 template<
typename TRAITS>                   
 
  493         const PHX::FieldTag & output,
 
  497   : use_descriptors_(true)
 
  506   accelerate_jacobian = 
false; 
 
  509   if(basis_dimension==2) {
 
  513   else if(basis_dimension==3) {
 
  521   this->addDependentField(dof_value);
 
  528 template<
typename TRAITS>                   
 
  533   this->utils.setFieldData(dof_value,fm);
 
  543 template<
typename TRAITS>                   
 
  550   if(!accelerate_jacobian) {
 
  552       EvaluateCurlWithSens_Vector<ScalarT,typename BasisValues2<double>::Array_CellBasisIPDim,3> functor(dof_value,
dof_curl_vector,basisValues.
curl_basis_vector);
 
  553       Kokkos::parallel_for(workset.num_cells,functor);
 
  557       Kokkos::parallel_for(workset.num_cells,functor);
 
  565       EvaluateCurlFastSens_Vector<ScalarT,typename BasisValues2<double>::Array_CellBasisIPDim,3> functor(dof_value,
dof_curl_vector,offsets_array,basisValues.
curl_basis_vector);
 
  566       Kokkos::parallel_for(workset.num_cells,functor);
 
  569       EvaluateCurlFastSens_Scalar<ScalarT,typename BasisValues2<double>::Array_CellBasisIP> functor(dof_value,
dof_curl_scalar,offsets_array,basisValues.
curl_basis_scalar);
 
  570       Kokkos::parallel_for(workset.num_cells,functor);
 
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
 
panzer::BasisDescriptor bd_
 
T & get(const std::string &name, T def_value)
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
Array_CellBasisIPDim curl_basis_vector
 
PHX::MDField< ScalarT, Cell, Point > dof_curl_scalar
 
PHX::MDField< ScalarT, Cell, Point, Dim > dof_curl_vector
 
const std::string & getType() const 
Get type of basis. 
 
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. 
 
DOFCurl(const Teuchos::ParameterList &p)
 
panzer::IntegrationDescriptor id_
 
bool isType(const std::string &name) const 
 
WorksetDetailsAccessor wda
 
#define TEUCHOS_ASSERT(assertion_test)
 
void evaluateFields(typename TRAITS::EvalData d)
 
PHX::MDField< const ScalarT, Cell, Point > dof_value
 
PHX::MDField< ScalarT, Cell, Point, Dim > dof_curl
 
Kokkos::View< const int *, PHX::Device > offsets
 
Array_CellBasisIP curl_basis_scalar