11#ifndef PANZER_DOF_CURL_IMPL_HPP
12#define PANZER_DOF_CURL_IMPL_HPP
17#include "Intrepid2_FunctionSpaceTools.hpp"
18#include "Phalanx_KokkosDeviceTypes.hpp"
25template <
typename ScalarT,
typename Array,
int spaceDim>
26class EvaluateCurlWithSens_Vector {
35 typedef typename PHX::Device execution_space;
37 EvaluateCurlWithSens_Vector(PHX::MDField<const ScalarT,Cell,Point> in_dof_value,
38 PHX::MDField<ScalarT,Cell,Point,Dim> in_dof_curl,
45 KOKKOS_INLINE_FUNCTION
46 void operator()(
const unsigned int cell)
const
49 for (
int d=0; d<spaceDim; d++) {
60template <
typename ScalarT,
typename ArrayT>
61void evaluateCurl_withSens_vector(
int numCells,
62 PHX::MDField<ScalarT,Cell,Point,Dim> &
dof_curl,
63 PHX::MDField<const ScalarT,Cell,Point> &
dof_value,
72 for (
int cell=0; cell<numCells; cell++) {
74 for (
int d=0; d<spaceDim; d++) {
87template <
typename ScalarT,
typename Array>
88class EvaluateCurlWithSens_Scalar {
89 PHX::MDField<const ScalarT,Cell,Point>
dof_value;
90 PHX::MDField<ScalarT,Cell,Point>
dof_curl;
97 typedef typename PHX::Device execution_space;
99 EvaluateCurlWithSens_Scalar(PHX::MDField<const ScalarT,Cell,Point> in_dof_value,
100 PHX::MDField<ScalarT,Cell,Point> in_dof_curl,
107 KOKKOS_INLINE_FUNCTION
108 void operator()(
const unsigned int cell)
const
120template <
typename ScalarT,
typename ArrayT>
121void evaluateCurl_withSens_scalar(
int numCells,
122 PHX::MDField<ScalarT,Cell,Point> &
dof_curl,
123 PHX::MDField<const ScalarT,Cell,Point> &
dof_value,
131 for (
int cell=0; cell<numCells; cell++) {
144template <
typename ScalarT,
typename Array,
int spaceDim>
145class EvaluateCurlFastSens_Vector {
146 PHX::MDField<const ScalarT,Cell,Point>
dof_value;
147 PHX::MDField<ScalarT,Cell,Point,Dim>
dof_curl;
155 typedef typename PHX::Device execution_space;
157 EvaluateCurlFastSens_Vector(PHX::MDField<const ScalarT,Cell,Point> in_dof_value,
158 PHX::MDField<ScalarT,Cell,Point,Dim> in_dof_curl,
159 PHX::View<const int*> in_offsets,
166 KOKKOS_INLINE_FUNCTION
167 void operator()(
const unsigned int cell)
const
170 for (
int d=0; d<spaceDim; d++) {
183template <
typename ScalarT,
typename ArrayT>
184void evaluateCurl_fastSens_vector(
int numCells,
185 PHX::MDField<ScalarT,Cell,Point,Dim> &
dof_curl,
186 PHX::MDField<const ScalarT,Cell,Point> &
dof_value,
187 const std::vector<int> &
offsets,
195 for (
int cell=0; cell<numCells; cell++) {
197 for (
int d=0; d<spaceDim; d++) {
213template <
typename ScalarT,
typename Array>
214class EvaluateCurlFastSens_Scalar {
215 PHX::MDField<const ScalarT,Cell,Point>
dof_value;
216 PHX::MDField<ScalarT,Cell,Point>
dof_curl;
224 typedef typename PHX::Device execution_space;
226 EvaluateCurlFastSens_Scalar(PHX::MDField<const ScalarT,Cell,Point> in_dof_value,
227 PHX::MDField<ScalarT,Cell,Point> in_dof_curl,
228 PHX::View<const int*> in_offsets,
235 KOKKOS_INLINE_FUNCTION
236 void operator()(
const unsigned int cell)
const
250template <
typename ScalarT,
typename ArrayT>
251void evaluateCurl_fastSens_scalar(
int numCells,
252 PHX::MDField<ScalarT,Cell,Point> &
dof_curl,
253 PHX::MDField<const ScalarT,Cell,Point> &
dof_value,
254 const std::vector<int> &
offsets,
261 for (
int cell=0; cell<numCells; cell++) {
285template<
typename EvalT,
typename TRAITS>
287DOFCurl(
const Teuchos::ParameterList & p) :
288 use_descriptors_(false),
293 Teuchos::RCP<const PureBasis> basis
294 = p.get< Teuchos::RCP<BasisIRLayout> >(
"Basis")->getBasis();
297 TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsCurl(),std::logic_error,
298 "DOFCurl: Basis of type \"" << basis->name() <<
"\" does not support CURL");
299 TEUCHOS_TEST_FOR_EXCEPTION(!basis->requiresOrientations(),std::logic_error,
300 "DOFCurl: Basis of type \"" << basis->name() <<
"\" in DOF Curl should require orientations. So we are throwing.");
305 dof_curl_scalar = PHX::MDField<ScalarT,Cell,Point>(p.get<std::string>(
"Curl Name"),
306 p.get< Teuchos::RCP<panzer::IntegrationRule> >(
"IR")->dl_scalar );
310 dof_curl_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(p.get<std::string>(
"Curl Name"),
311 p.get< Teuchos::RCP<panzer::IntegrationRule> >(
"IR")->dl_vector );
315 { TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
"DOFCurl only works for 2D and 3D basis functions"); }
325template<
typename EvalT,
typename TRAITS>
327DOFCurl(
const PHX::FieldTag & input,
328 const PHX::FieldTag & output,
332 : use_descriptors_(true)
351 { TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
"DOFCurl only works for 2D and 3D basis functions"); }
361template<
typename EvalT,
typename TRAITS>
367 if(basis_dimension==3)
368 this->utils.setFieldData(dof_curl_vector,fm);
370 this->utils.setFieldData(dof_curl_scalar,fm);
372 if(not use_descriptors_)
377template<
typename EvalT,
typename TRAITS>
382 : *this->wda(workset).bases[basis_index];
384 if(basis_dimension==3) {
385 EvaluateCurlWithSens_Vector<ScalarT,typename BasisValues2<double>::Array_CellBasisIPDim,3> functor(
dof_value,dof_curl_vector,basisValues.
curl_basis_vector);
386 Kokkos::parallel_for(workset.num_cells,functor);
389 EvaluateCurlWithSens_Scalar<ScalarT,typename BasisValues2<double>::Array_CellBasisIP> functor(
dof_value,dof_curl_scalar,basisValues.
curl_basis_scalar);
390 Kokkos::parallel_for(workset.num_cells,functor);
401template<
typename TRAITS>
403DOFCurl(
const Teuchos::ParameterList & p) :
404 use_descriptors_(false),
409 Teuchos::RCP<const PureBasis> basis
410 = p.get< Teuchos::RCP<BasisIRLayout> >(
"Basis")->getBasis();
414 if(p.isType<Teuchos::RCP<
const std::vector<int> > >(
"Jacobian Offsets Vector")) {
415 offsets = *p.get<Teuchos::RCP<const std::vector<int> > >(
"Jacobian Offsets Vector");
418 PHX::View<int*> offsets_array_nc(
"offsets",
offsets.size());
419 auto offsets_array_nc_h = Kokkos::create_mirror_view(offsets_array_nc);
420 for(std::size_t i=0;i<
offsets.size();i++)
421 offsets_array_nc_h(i) =
offsets[i];
422 Kokkos::deep_copy(offsets_array_nc, offsets_array_nc_h);
423 offsets_array = offsets_array_nc;
425 accelerate_jacobian =
true;
428 accelerate_jacobian =
false;
431 TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsCurl(),std::logic_error,
432 "DOFCurl: Basis of type \"" << basis->name() <<
"\" does not support CURL");
433 TEUCHOS_TEST_FOR_EXCEPTION(!basis->requiresOrientations(),std::logic_error,
434 "DOFCurl: Basis of type \"" << basis->name() <<
"\" in DOF Curl should require orientations. So we are throwing.");
439 dof_curl_scalar = PHX::MDField<ScalarT,Cell,Point>(p.get<std::string>(
"Curl Name"),
440 p.get< Teuchos::RCP<panzer::IntegrationRule> >(
"IR")->dl_scalar );
444 dof_curl_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(p.get<std::string>(
"Curl Name"),
445 p.get< Teuchos::RCP<panzer::IntegrationRule> >(
"IR")->dl_vector );
449 { TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
"DOFCurl only works for 2D and 3D basis functions"); }
459template<
typename TRAITS>
461DOFCurl(
const PHX::FieldTag & input,
462 const PHX::FieldTag & output,
466 : use_descriptors_(true)
475 accelerate_jacobian =
false;
487 { TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
"DOFCurl only works for 2D and 3D basis functions"); }
497template<
typename TRAITS>
503 if(basis_dimension==3)
504 this->utils.setFieldData(dof_curl_vector,fm);
506 this->utils.setFieldData(dof_curl_scalar,fm);
508 if(not use_descriptors_)
512template<
typename TRAITS>
517 : *this->wda(workset).bases[basis_index];
519 if(!accelerate_jacobian) {
520 if(basis_dimension==3) {
523 EvaluateCurlWithSens_Vector<ScalarT,Array,3> functor(
dof_value,dof_curl_vector,curl_basis_vector);
524 Kokkos::parallel_for(workset.num_cells,functor);
529 EvaluateCurlWithSens_Scalar<ScalarT,Array> functor(
dof_value,dof_curl_scalar,curl_basis_scalar);
530 Kokkos::parallel_for(workset.num_cells,functor);
537 if(basis_dimension==3) {
540 EvaluateCurlFastSens_Vector<ScalarT,Array,3> functor(
dof_value,dof_curl_vector,offsets_array,curl_basis_vector);
541 Kokkos::parallel_for(workset.num_cells,functor);
546 EvaluateCurlFastSens_Scalar<ScalarT,Array> functor(
dof_value,dof_curl_scalar,offsets_array,curl_basis_scalar);
547 Kokkos::parallel_for(workset.num_cells,functor);
PHX::MDField< ScalarT, Cell, Point, Dim > dof_curl
PHX::View< const int * > offsets
PHX::MDField< const ScalarT, Cell, Point > dof_value
const std::string & getType() const
Get type of basis.
PHX::MDField< const Scalar, Cell, BASIS, IP > ConstArray_CellBasisIP
Array_CellBasisIPDim curl_basis_vector
ConstArray_CellBasisIPDim getCurlVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the curl of a 3D vector basis evaluated at mesh points.
Array_CellBasisIP curl_basis_scalar
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
PHX::MDField< const Scalar, Cell, BASIS, IP, Dim > ConstArray_CellBasisIPDim
ConstArray_CellBasisIP getCurl2DVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the curl of a 2D vector basis evaluated at mesh points.
DOFCurl(const Teuchos::ParameterList &p)
PHX::MDField< ScalarT, Cell, Point > dof_curl_scalar
void evaluateFields(typename TRAITS::EvalData d)
PHX::MDField< const ScalarT, Cell, Point > dof_value
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
panzer::BasisDescriptor bd_
PHX::MDField< ScalarT, Cell, Point, Dim > dof_curl_vector
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.