77 KOKKOS_INLINE_FUNCTION
78 void operator()(
const unsigned int cell)
const
80 for (
int pt=0; pt<numPoints; pt++) {
81 for (
int d=0; d<spaceDim; d++) {
84 dof_curl(cell,pt,d) =
dof_value(cell, 0) * curl_basis(cell, 0, pt, d);
85 for (
int bf=1; bf<numFields; bf++)
86 dof_curl(cell,pt,d) +=
dof_value(cell, bf) * curl_basis(cell, bf, pt, d);
92template <
typename ScalarT,
typename ArrayT>
93void evaluateCurl_withSens_vector(
int numCells,
94 PHX::MDField<ScalarT,Cell,Point,Dim> & dof_curl,
95 PHX::MDField<const ScalarT,Cell,Point> & dof_value,
96 const ArrayT & curl_basis)
100 int numFields = curl_basis.extent(1);
101 int numPoints = curl_basis.extent(2);
102 int spaceDim = curl_basis.extent(3);
104 for (
int cell=0; cell<numCells; cell++) {
105 for (
int pt=0; pt<numPoints; pt++) {
106 for (
int d=0; d<spaceDim; d++) {
109 dof_curl(cell,pt,d) = dof_value(cell, 0) * curl_basis(cell, 0, pt, d);
110 for (
int bf=1; bf<numFields; bf++)
111 dof_curl(cell,pt,d) += dof_value(cell, bf) * curl_basis(cell, bf, pt, d);
119template <
typename ScalarT,
typename Array>
120class 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
142 for (
int pt=0; pt<numPoints; pt++) {
145 dof_curl(cell,pt) = dof_value(cell, 0) * curl_basis(cell, 0, pt);
146 for (
int bf=1; bf<numFields; bf++)
147 dof_curl(cell,pt) += dof_value(cell, bf) * curl_basis(cell, bf, pt);
152template <
typename ScalarT,
typename ArrayT>
153void 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++) {
164 for (
int pt=0; pt<numPoints; pt++) {
167 dof_curl(cell,pt) = dof_value(cell, 0) * curl_basis(cell, 0, pt);
168 for (
int bf=1; bf<numFields; bf++)
169 dof_curl(cell,pt) += dof_value(cell, bf) * curl_basis(cell, bf, pt);
176template <
typename ScalarT,
typename Array,
int spaceDim>
177class EvaluateCurlFastSens_Vector {
178 PHX::MDField<const ScalarT,Cell,Point> dof_value;
179 PHX::MDField<ScalarT,Cell,Point,Dim> dof_curl;
180 PHX::View<const int*> offsets;
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 PHX::View<const int*> 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
201 for (
int pt=0; pt<numPoints; pt++) {
202 for (
int d=0; d<spaceDim; d++) {
205 dof_curl(cell,pt,d) = dof_value(cell, 0).val() * curl_basis(cell, 0, pt, d);
206 dof_curl(cell,pt,d).fastAccessDx(offsets(0)) = dof_value(cell, 0).fastAccessDx(offsets(0)) * curl_basis(cell, 0, pt, d);
207 for (
int bf=1; bf<numFields; bf++) {
208 dof_curl(cell,pt,d).val() += dof_value(cell, bf).val() * curl_basis(cell, bf, pt, d);
209 dof_curl(cell,pt,d).fastAccessDx(offsets(bf)) += dof_value(cell, bf).fastAccessDx(offsets(bf)) * curl_basis(cell, bf, pt, d);
215template <
typename ScalarT,
typename ArrayT>
216void 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++) {
228 for (
int pt=0; pt<numPoints; pt++) {
229 for (
int d=0; d<spaceDim; d++) {
232 dof_curl(cell,pt,d) = ScalarT(numFields, dof_value(cell, 0).val() * curl_basis(cell, 0, pt, d));
233 dof_curl(cell,pt,d).fastAccessDx(offsets[0]) = dof_value(cell, 0).fastAccessDx(offsets[0]) * curl_basis(cell, 0, pt, d);
234 for (
int bf=1; bf<numFields; bf++) {
235 dof_curl(cell,pt,d).val() += dof_value(cell, bf).val() * curl_basis(cell, bf, pt, d);
236 dof_curl(cell,pt,d).fastAccessDx(offsets[bf]) += dof_value(cell, bf).fastAccessDx(offsets[bf]) * curl_basis(cell, bf, pt, d);
245template <
typename ScalarT,
typename Array>
246class EvaluateCurlFastSens_Scalar {
247 PHX::MDField<const ScalarT,Cell,Point> dof_value;
248 PHX::MDField<ScalarT,Cell,Point> dof_curl;
249 PHX::View<const int*> 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 PHX::View<const int*> 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
270 for (
int pt=0; pt<numPoints; pt++) {
273 dof_curl(cell,pt) = dof_value(cell, 0).val() * curl_basis(cell, 0, pt);
274 dof_curl(cell,pt).fastAccessDx(offsets(0)) = dof_value(cell, 0).fastAccessDx(offsets(0)) * curl_basis(cell, 0, pt);
275 for (
int bf=1; bf<numFields; bf++) {
276 dof_curl(cell,pt).val() += dof_value(cell, bf).val() * curl_basis(cell, bf, pt);
277 dof_curl(cell,pt).fastAccessDx(offsets(bf)) += dof_value(cell, bf).fastAccessDx(offsets(bf)) * curl_basis(cell, bf, pt);
282template <
typename ScalarT,
typename ArrayT>
283void 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++) {
294 for (
int pt=0; pt<numPoints; pt++) {
297 dof_curl(cell,pt) = ScalarT(numFields, dof_value(cell, 0).val() * curl_basis(cell, 0, pt));
298 dof_curl(cell,pt).fastAccessDx(offsets[0]) = dof_value(cell, 0).fastAccessDx(offsets[0]) * curl_basis(cell, 0, pt);
299 for (
int bf=1; bf<numFields; bf++) {
300 dof_curl(cell,pt).val() += dof_value(cell, bf).val() * curl_basis(cell, bf, pt);
301 dof_curl(cell,pt).fastAccessDx(offsets[bf]) += dof_value(cell, bf).fastAccessDx(offsets[bf]) * curl_basis(cell, bf, pt);
317template<
typename EvalT,
typename TRAITS>
319DOFCurl(
const Teuchos::ParameterList & p) :
325 Teuchos::RCP<const PureBasis> basis
326 = p.get< Teuchos::RCP<BasisIRLayout> >(
"Basis")->getBasis();
329 TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsCurl(),std::logic_error,
330 "DOFCurl: Basis of type \"" << basis->name() <<
"\" does not support CURL");
331 TEUCHOS_TEST_FOR_EXCEPTION(!basis->requiresOrientations(),std::logic_error,
332 "DOFCurl: Basis of type \"" << basis->name() <<
"\" in DOF Curl should require orientations. So we are throwing.");
337 dof_curl_scalar = PHX::MDField<ScalarT,Cell,Point>(p.get<std::string>(
"Curl Name"),
338 p.get< Teuchos::RCP<panzer::IntegrationRule> >(
"IR")->dl_scalar );
342 dof_curl_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(p.get<std::string>(
"Curl Name"),
343 p.get< Teuchos::RCP<panzer::IntegrationRule> >(
"IR")->dl_vector );
347 { TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
"DOFCurl only works for 2D and 3D basis functions"); }
357template<
typename EvalT,
typename TRAITS>
359DOFCurl(
const PHX::FieldTag & input,
360 const PHX::FieldTag & output,
369 TEUCHOS_ASSERT(
bd_.getType()==
"HCurl");
383 { TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
"DOFCurl only works for 2D and 3D basis functions"); }
393template<
typename EvalT,
typename TRAITS>
409template<
typename EvalT,
typename TRAITS>
418 Kokkos::parallel_for(workset.num_cells,functor);
422 Kokkos::parallel_for(workset.num_cells,functor);