Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_BasisValues2_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #include "PanzerDiscFE_config.hpp"
44 #include "Panzer_Traits.hpp"
45 
47 #include "Kokkos_ViewFactory.hpp"
48 
49 #include "Intrepid2_Utils.hpp"
50 #include "Intrepid2_FunctionSpaceTools.hpp"
51 #include "Intrepid2_Orientation.hpp"
52 #include "Intrepid2_OrientationTools.hpp"
53 
54 
55 namespace panzer {
56 
57 template <typename Scalar>
59 evaluateValues(const PHX::MDField<Scalar,IP,Dim,void,void,void,void,void,void> & cub_points,
60  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
61  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
62  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
63  const int in_num_cells)
64 {
65  PHX::MDField<Scalar,Cell,IP> weighted_measure;
66  PHX::MDField<Scalar,Cell,NODE,Dim> vertex_coordinates;
67  build_weighted = false;
68  evaluateValues(cub_points,jac,jac_det,jac_inv,weighted_measure,vertex_coordinates,false,in_num_cells);
69 }
70 
71 template <typename Scalar>
73 evaluateValues(const PHX::MDField<Scalar,IP,Dim,void,void,void,void,void,void> & cub_points,
74  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
75  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
76  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
77  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
78  const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
79  bool use_vertex_coordinates,
80  const int in_num_cells)
81 {
82  MDFieldArrayFactory af("",ddims_,true);
83 
84  int num_dim = basis_layout->dimension();
85 
86  // currently this just copies on the basis objects are converted
87  // in intrepid
88  evaluateReferenceValues(cub_points,compute_derivatives,use_vertex_coordinates);
89 
90  const int num_cells = in_num_cells < 0 ? jac.extent(0) : in_num_cells;
91  const std::pair<int,int> cell_range(0,num_cells);
92 
93  PureBasis::EElementSpace elmtspace = getElementSpace();
94  if(elmtspace==PureBasis::CONST ||
95  elmtspace==PureBasis::HGRAD) {
96  auto s_basis_scalar = Kokkos::subview(basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
97  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
98  HGRADtransformVALUE(s_basis_scalar,
99  basis_ref_scalar.get_view());
100 
101  if(build_weighted) {
102  auto s_weighted_basis_scalar = Kokkos::subview(weighted_basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
103  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
104  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
105  multiplyMeasure(s_weighted_basis_scalar,
106  s_weighted_measure,
107  s_basis_scalar);
108  }
109  }
110  else if(elmtspace==PureBasis::HCURL) {
111  auto s_basis_vector = Kokkos::subview(basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
112  auto s_jac_inv = Kokkos::subview(jac_inv.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
113  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
114  HCURLtransformVALUE(s_basis_vector,
115  s_jac_inv,
116  basis_ref_vector.get_view());
117 
118  if(build_weighted) {
119  auto s_weighted_basis_vector = Kokkos::subview(weighted_basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
120  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
121  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
122  multiplyMeasure(s_weighted_basis_vector,
123  s_weighted_measure,
124  s_basis_vector);
125  }
126  }
127  else if(elmtspace==PureBasis::HDIV)
128  {
129  auto s_basis_vector = Kokkos::subview(basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
130  auto s_jac = Kokkos::subview(jac.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
131  auto s_jac_det = Kokkos::subview(jac_det.get_view(),cell_range,Kokkos::ALL());
132  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
133  HDIVtransformVALUE(s_basis_vector,
134  s_jac,
135  s_jac_det,
136  basis_ref_vector.get_view());
137 
138  if(build_weighted) {
139  auto s_weighted_basis_vector = Kokkos::subview(weighted_basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
140  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
141  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
142  multiplyMeasure(s_weighted_basis_vector,
143  s_weighted_measure,
144  s_basis_vector);
145  }
146  }
147  else if(elmtspace==PureBasis::HVOL)
148  {
149  auto s_basis_scalar = Kokkos::subview(basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
150  auto s_jac_det = Kokkos::subview(jac_det.get_view(),cell_range,Kokkos::ALL());
151  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
152  HVOLtransformVALUE(s_basis_scalar,
153  s_jac_det,
154  basis_ref_scalar.get_view());
155 
156  if(build_weighted) {
157  auto s_weighted_basis_scalar = Kokkos::subview(weighted_basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
158  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
159  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
160  multiplyMeasure(s_weighted_basis_scalar,
161  s_weighted_measure,
162  s_basis_scalar);
163  }
164  }
165  else { TEUCHOS_ASSERT(false); }
166 
167  if(elmtspace==PureBasis::HGRAD && compute_derivatives) {
168  auto s_grad_basis = Kokkos::subview(grad_basis.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
169  auto s_jac_inv = Kokkos::subview(jac_inv.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
170  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
171  HGRADtransformGRAD(s_grad_basis,
172  s_jac_inv,
173  grad_basis_ref.get_view());
174 
175  if(build_weighted) {
176  auto s_weighted_grad_basis = Kokkos::subview(weighted_grad_basis.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
177  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
178  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
179  multiplyMeasure(s_weighted_grad_basis,
180  s_weighted_measure,
181  s_grad_basis);
182  }
183  }
184  else if(elmtspace==PureBasis::HCURL && num_dim==2 && compute_derivatives) {
185  auto s_curl_basis_scalar = Kokkos::subview(curl_basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
186  auto s_jac_det = Kokkos::subview(jac_det.get_view(),cell_range,Kokkos::ALL());
187  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
188  HDIVtransformDIV(s_curl_basis_scalar,
189  s_jac_det, // note only volume deformation is needed!
190  // this relates directly to this being in
191  // the divergence space in 2D!
192  curl_basis_ref_scalar.get_view());
193 
194  if(build_weighted) {
195  auto s_weighted_curl_basis_scalar = Kokkos::subview(weighted_curl_basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
196  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
197  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
198  multiplyMeasure(s_weighted_curl_basis_scalar,
199  s_weighted_measure,
200  s_curl_basis_scalar);
201  }
202  }
203  else if(elmtspace==PureBasis::HCURL && num_dim==3 && compute_derivatives) {
204  auto s_curl_basis_vector = Kokkos::subview(curl_basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
205  auto s_jac = Kokkos::subview(jac.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
206  auto s_jac_det = Kokkos::subview(jac_det.get_view(),cell_range,Kokkos::ALL());
207  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
208  HCURLtransformCURL(s_curl_basis_vector,
209  s_jac,
210  s_jac_det,
211  curl_basis_ref_vector.get_view());
212 
213  if(build_weighted) {
214  auto s_weighted_curl_basis_vector = Kokkos::subview(weighted_curl_basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
215  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
216  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
217  multiplyMeasure(s_weighted_curl_basis_vector,
218  s_weighted_measure,
219  s_curl_basis_vector);
220  }
221  }
222  else if(elmtspace==PureBasis::HDIV && compute_derivatives) {
223  auto s_div_basis = Kokkos::subview(div_basis.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
224  auto s_jac_det = Kokkos::subview(jac_det.get_view(),cell_range,Kokkos::ALL());
225  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
226  HDIVtransformDIV(s_div_basis,
227  s_jac_det,
228  div_basis_ref.get_view());
229 
230  if(build_weighted) {
231  auto s_weighted_div_basis = Kokkos::subview(weighted_div_basis.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
232  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
233  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::
234  multiplyMeasure(s_weighted_div_basis,
235  s_weighted_measure,
236  s_div_basis);
237  }
238  }
239 
240  // If basis supports coordinate values at basis points, then
241  // compute these values
242  if(use_vertex_coordinates) {
243  // Teuchos::RCP<Intrepid2::DofCoordsInterface<ArrayDynamic> > coords
244  // = Teuchos::rcp_dynamic_cast<Intrepid2::DofCoordsInterface<ArrayDynamic> >(intrepid_basis);
245  // if (!Teuchos::is_null(coords)) {
246 /*
247  ArrayDynamic dyn_basis_coordinates_ref = af.buildArray<Scalar,BASIS,Dim>("basis_coordinates_ref",basis_coordinates_ref.extent(0),basis_coordinates_ref.extent(1));
248  coords->getDofCoords(dyn_basis_coordinates_ref);
249 
250  // fill in basis coordinates
251  for (size_type i = 0; i < basis_coordinates_ref.extent(0); ++i)
252  for (size_type j = 0; j < basis_coordinates_ref.extent(1); ++j)
253  basis_coordinates_ref(i,j) = dyn_basis_coordinates_ref(i,j);
254 */
255 
256  auto s_basis_coordinates = Kokkos::subview(basis_coordinates.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
257  auto s_vertex_coordinates = Kokkos::subview(vertex_coordinates.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
258  Intrepid2::CellTools<PHX::Device::execution_space> cell_tools;
259  cell_tools.mapToPhysicalFrame(s_basis_coordinates,
260  basis_coordinates_ref.get_view(),
261  s_vertex_coordinates,
262  intrepid_basis->getBaseCellTopology());
263  }
264 }
265 
266 
267 
268 
269 
270 
271 
272 template <typename Scalar>
274 evaluateBasisCoordinates(const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
275  const int in_num_cells)
276 {
277  MDFieldArrayFactory af("",ddims_,true);
278 
279  // intrepid_basis->getDofCoords requires DynRankView, but basis_coordinates_ref is more of a static View
280  // We use an auxiliary 'dyn' array to get around this
281  using coordsScalarType = typename Intrepid2::Basis<PHX::Device::execution_space,Scalar,Scalar>::scalarType;
282  auto dyn_basis_coordinates_ref = af.buildArray<coordsScalarType,BASIS,Dim>("basis_coordinates_ref",
283  basis_coordinates_ref.extent(0),
284  basis_coordinates_ref.extent(1));
285  intrepid_basis->getDofCoords(dyn_basis_coordinates_ref.get_view());
286 
287  // fill in basis coordinates
288  for (int i = 0; i < basis_coordinates_ref.extent_int(0); ++i)
289  for (int j = 0; j < basis_coordinates_ref.extent_int(1); ++j)
290  basis_coordinates_ref(i,j) = dyn_basis_coordinates_ref(i,j);
291 
292  const int num_cells = in_num_cells < 0 ? vertex_coordinates.extent(0) : in_num_cells;
293  const std::pair<int,int> cell_range(0,num_cells);
294  const auto s_basis_coordinates = Kokkos::subview(basis_coordinates.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
295  const auto s_vertex_coordinates = Kokkos::subview(vertex_coordinates.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
296 
297  Intrepid2::CellTools<PHX::Device::execution_space> cell_tools;
298  cell_tools.mapToPhysicalFrame(s_basis_coordinates,
299  basis_coordinates_ref.get_view(),
300  s_vertex_coordinates,
301  intrepid_basis->getBaseCellTopology());
302 }
303 
304 
305 
306 
307 
308 template <typename Scalar>
310 evaluateValues(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
311  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
312  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
313  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
314  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
315  const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
316  bool use_vertex_coordinates,
317  const int in_num_cells)
318 {
319 
320  PureBasis::EElementSpace elmtspace = getElementSpace();
321 
322  if(elmtspace == PureBasis::CONST){
323  evaluateValues_Const(cub_points,jac_inv,weighted_measure,in_num_cells);
324  } else if(elmtspace == PureBasis::HVOL){
325  evaluateValues_HVol(cub_points,jac_det,jac_inv,weighted_measure,in_num_cells);
326  } else if(elmtspace == PureBasis::HGRAD){
327  evaluateValues_HGrad(cub_points,jac_inv,weighted_measure,in_num_cells);
328  } else if(elmtspace == PureBasis::HCURL){
329  evaluateValues_HCurl(cub_points,jac,jac_det,jac_inv,weighted_measure,in_num_cells);
330  } else if(elmtspace == PureBasis::HDIV){
331  evaluateValues_HDiv(cub_points,jac,jac_det,weighted_measure,in_num_cells);
332  } else {
333  TEUCHOS_TEST_FOR_EXCEPT_MSG(true,"panzer::BasisValues2::evaluateValues : Element space not recognized.");
334  }
335 
336  if(use_vertex_coordinates) {
337  TEUCHOS_TEST_FOR_EXCEPT_MSG(elmtspace == PureBasis::CONST,"panzer::BasisValues2::evaluateValues : Const basis cannot have basis coordinates.");
338  evaluateBasisCoordinates(vertex_coordinates);
339  }
340 
341 }
342 
343 template <typename Scalar>
345 evaluateValues_Const(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
346  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
347  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
348  const int in_num_cells)
349 {
350 
351  TEUCHOS_ASSERT(getElementSpace() == PureBasis::CONST);
352 
353  typedef Intrepid2::FunctionSpaceTools<PHX::Device::execution_space> fst;
354  MDFieldArrayFactory af("",ddims_,true);
355 
356  const panzer::PureBasis & basis = *(basis_layout->getBasis());
357 
358  const int num_points = basis_layout->numPoints();
359  const int num_basis = basis.cardinality();
360  const int num_dim = basis_layout->dimension();
361  const int num_cells = in_num_cells < 0 ? basis_layout->numCells() : in_num_cells;
362 
363  auto cell_basis_scalar = af.buildStaticArray<Scalar,Cell,BASIS,IP>("cell_basis_scalar",1,num_basis,num_points);
364  auto cell_cub_points = af.buildStaticArray<Scalar,IP,Dim>("cell_cub_points",num_points,num_dim);
365  auto cell_grad_basis = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("cell_grad_basis",1,num_basis,num_points,num_dim);
366  auto cell_jac_inv = af.buildStaticArray<Scalar,Cell,IP,Dim,Dim>("cell_jac_inv",1,num_points,num_dim,num_dim);
367 
368  auto cell_basis_ref_scalar = af.buildStaticArray<Scalar,BASIS,IP>("cell_basis_ref_scalar",num_basis,num_points);
369  auto cell_grad_basis_ref = af.buildStaticArray<Scalar,BASIS,IP,Dim>("cell_grad_basis_ref",num_basis,num_points,num_dim);
370 
371  for(int cell=0;cell<num_cells;++cell){
372 
373  // =============================================
374  // Load external into cell-local arrays
375 
376  for(int p=0;p<num_points;++p)
377  for(int d=0;d<num_dim;++d)
378  for(int d2=0;d2<num_dim;++d2)
379  cell_jac_inv(0,p,d,d2)=jac_inv(cell,p,d,d2);
380  for(int p=0;p<num_points;++p)
381  for(int d=0;d<num_dim;++d)
382  cell_cub_points(p,d)=cub_points(cell,p,d);
383 
384  // =============================================
385  // Load Reference Values
386 
387  intrepid_basis->getValues(cell_basis_ref_scalar.get_view(),cell_cub_points.get_view(),Intrepid2::OPERATOR_VALUE);
388 
389  if(compute_derivatives){
390  Kokkos::deep_copy(cell_grad_basis_ref.get_view(),0.0);
391  }
392 
393  // =============================================
394  // Transform reference values to physical values
395 
396  fst::HGRADtransformVALUE(cell_basis_scalar.get_view(),cell_basis_ref_scalar.get_view());
397  for(int b=0;b<num_basis;++b)
398  for(int p=0;p<num_points;++p)
399  basis_scalar(cell,b,p)=cell_basis_scalar(0,b,p);
400 
401  if(compute_derivatives){
402  fst::HGRADtransformGRAD(cell_grad_basis.get_view(),cell_jac_inv.get_view(),cell_grad_basis_ref.get_view());
403  for(int b=0;b<num_basis;++b)
404  for(int p=0;p<num_points;++p)
405  for(int d=0;d<num_dim;++d)
406  grad_basis(cell,b,p,d)=cell_grad_basis(0,b,p,d);
407  }
408  // =============================================
409  }
410 
411 
412  if(build_weighted){
413  const std::pair<int,int> cell_range(0,num_cells);
414  auto s_weighted_basis_scalar = Kokkos::subview(weighted_basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
415  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
416  auto s_basis_scalar = Kokkos::subview(basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
417  fst::multiplyMeasure(s_weighted_basis_scalar,s_weighted_measure,s_basis_scalar);
418  if(compute_derivatives){
419  auto s_weighted_grad_basis = Kokkos::subview(weighted_grad_basis.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
420  auto s_grad_basis = Kokkos::subview(grad_basis.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
421  fst::multiplyMeasure(s_weighted_grad_basis,s_weighted_measure,s_grad_basis);
422  }
423  }
424 
425 
426 }
427 
428 template <typename Scalar>
430 evaluateValues_HVol(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
431  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
432  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & /* jac_inv */,
433  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
434  const int in_num_cells)
435 {
436 
437  TEUCHOS_ASSERT(getElementSpace() == PureBasis::HVOL);
438 
439  typedef Intrepid2::FunctionSpaceTools<PHX::Device::execution_space> fst;
440  MDFieldArrayFactory af("",ddims_,true);
441 
442  const panzer::PureBasis & basis = *(basis_layout->getBasis());
443 
444  const int num_points = basis_layout->numPoints();
445  const int num_basis = basis.cardinality();
446  const int num_dim = basis_layout->dimension();
447  const int num_cells = in_num_cells < 0 ? basis_layout->numCells() : in_num_cells;
448 
449  auto cell_basis_scalar = af.buildStaticArray<Scalar,Cell,BASIS,IP>("cell_basis_scalar",1,num_basis,num_points);
450  auto cell_cub_points = af.buildStaticArray<Scalar,IP,Dim>("cell_cub_points",num_points,num_dim);
451  auto cell_jac_det = af.buildStaticArray<Scalar,Cell,IP>("cell_jac_det",1,num_points);
452 
453  auto cell_basis_ref_scalar = af.buildStaticArray<Scalar,BASIS,IP>("cell_basis_ref_scalar",num_basis,num_points);
454 
455  for(int cell=0;cell<num_cells;++cell){
456 
457  // =============================================
458  // Load external into cell-local arrays
459 
460  for(int p=0;p<num_points;++p)
461  cell_jac_det(0,p)=jac_det(cell,p);
462  for(int p=0;p<num_points;++p)
463  for(int d=0;d<num_dim;++d)
464  cell_cub_points(p,d)=cub_points(cell,p,d);
465 
466  // =============================================
467  // Load Reference Values
468 
469  intrepid_basis->getValues(cell_basis_ref_scalar.get_view(),cell_cub_points.get_view(),Intrepid2::OPERATOR_VALUE);
470 
471  // =============================================
472  // Transform reference values to physical values
473 
474  fst::HVOLtransformVALUE(cell_basis_scalar.get_view(),cell_jac_det.get_view(),cell_basis_ref_scalar.get_view());
475  for(int b=0;b<num_basis;++b)
476  for(int p=0;p<num_points;++p)
477  basis_scalar(cell,b,p)=cell_basis_scalar(0,b,p);
478 
479  // =============================================
480  }
481 
482 
483  if(build_weighted) {
484  const std::pair<int,int> cell_range(0,num_cells);
485  auto s_weighted_basis_scalar = Kokkos::subview(weighted_basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
486  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
487  auto s_basis_scalar = Kokkos::subview(basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
488  fst::multiplyMeasure(s_weighted_basis_scalar,s_weighted_measure,s_basis_scalar);
489  }
490 
491 }
492 
493 template <typename Scalar>
495 evaluateValues_HGrad(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
496  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
497  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
498  const int in_num_cells)
499 {
500 
501  TEUCHOS_ASSERT(getElementSpace() == PureBasis::HGRAD);
502 
503  typedef Intrepid2::FunctionSpaceTools<PHX::Device::execution_space> fst;
504  MDFieldArrayFactory af("",ddims_,true);
505 
506  const panzer::PureBasis & basis = *(basis_layout->getBasis());
507 
508  const int num_points = basis_layout->numPoints();
509  const int num_basis = basis.cardinality();
510  const int num_dim = basis_layout->dimension();
511  const int num_cells = in_num_cells < 0 ? cub_points.extent(0) : in_num_cells;
512 
513  auto cell_basis_scalar = af.buildStaticArray<Scalar,Cell,BASIS,IP>("cell_basis_scalar",1,num_basis,num_points);
514  auto cell_cub_points = af.buildStaticArray<Scalar,IP,Dim>("cell_cub_points",num_points,num_dim);
515  auto cell_grad_basis = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("cell_grad_basis",1,num_basis,num_points,num_dim);
516  auto cell_jac_inv = af.buildStaticArray<Scalar,Cell,IP,Dim,Dim>("cell_jac_inv",1,num_points,num_dim,num_dim);
517 
518  auto cell_basis_ref_scalar = af.buildStaticArray<Scalar,BASIS,IP>("cell_basis_ref_scalar",num_basis,num_points);
519  auto cell_grad_basis_ref = af.buildStaticArray<Scalar,BASIS,IP,Dim>("cell_grad_basis_ref",num_basis,num_points,num_dim);
520 
521  for(int cell=0;cell<num_cells;++cell){
522 
523  // =============================================
524  // Load external into cell-local arrays
525 
526  for(int p=0;p<num_points;++p)
527  for(int d=0;d<num_dim;++d)
528  for(int d2=0;d2<num_dim;++d2)
529  cell_jac_inv(0,p,d,d2)=jac_inv(cell,p,d,d2);
530  for(int p=0;p<num_points;++p)
531  for(int d=0;d<num_dim;++d)
532  cell_cub_points(p,d)=cub_points(cell,p,d);
533 
534  // =============================================
535  // Load Reference Values
536 
537  intrepid_basis->getValues(cell_basis_ref_scalar.get_view(),cell_cub_points.get_view(),Intrepid2::OPERATOR_VALUE);
538 
539  if(compute_derivatives){
540  intrepid_basis->getValues(cell_grad_basis_ref.get_view(),cell_cub_points.get_view(),Intrepid2::OPERATOR_GRAD);
541  }
542 
543  // =============================================
544  // Transform reference values to physical values
545 
546  fst::HGRADtransformVALUE(cell_basis_scalar.get_view(),cell_basis_ref_scalar.get_view());
547  for(int b=0;b<num_basis;++b)
548  for(int p=0;p<num_points;++p)
549  basis_scalar(cell,b,p)=cell_basis_scalar(0,b,p);
550 
551  if(compute_derivatives){
552  fst::HGRADtransformGRAD(cell_grad_basis.get_view(),cell_jac_inv.get_view(),cell_grad_basis_ref.get_view());
553  for(int b=0;b<num_basis;++b)
554  for(int p=0;p<num_points;++p)
555  for(int d=0;d<num_dim;++d)
556  grad_basis(cell,b,p,d)=cell_grad_basis(0,b,p,d);
557  }
558  // =============================================
559  }
560 
561  if(build_weighted){
562  const std::pair<int,int> cell_range(0,num_cells);
563  auto s_weighted_basis_scalar = Kokkos::subview(weighted_basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
564  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
565  auto s_basis_scalar = Kokkos::subview(basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
566  fst::multiplyMeasure(s_weighted_basis_scalar,s_weighted_measure,s_basis_scalar);
567  if(compute_derivatives){
568  auto s_weighted_grad_basis = Kokkos::subview(weighted_grad_basis.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
569  auto s_grad_basis = Kokkos::subview(grad_basis.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
570  fst::multiplyMeasure(s_weighted_grad_basis,s_weighted_measure,s_grad_basis);
571  }
572  }
573 
574 }
575 
576 
577 template <typename Scalar>
579 evaluateValues_HCurl(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
580  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
581  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
582  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
583  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
584  const int in_num_cells)
585 {
586 
587  TEUCHOS_ASSERT(getElementSpace() == PureBasis::HCURL);
588 
589 
590  typedef Intrepid2::FunctionSpaceTools<PHX::Device::execution_space> fst;
591  MDFieldArrayFactory af("",ddims_,true);
592 
593  const panzer::PureBasis & basis = *(basis_layout->getBasis());
594 
595  const int num_points = basis_layout->numPoints();
596  const int num_basis = basis.cardinality();
597  const int num_dim = basis_layout->dimension();
598  const int num_cells = in_num_cells < 0 ? basis_layout->numCells() : in_num_cells;
599 
600  auto cell_cub_points = af.buildStaticArray<Scalar,IP,Dim>("cell_cub_points",num_points,num_dim);
601  auto cell_jac = af.buildStaticArray<Scalar,Cell,IP,Dim,Dim>("cell_jac",1,num_points,num_dim,num_dim);
602  auto cell_jac_inv = af.buildStaticArray<Scalar,Cell,IP,Dim,Dim>("cell_jac_inv",1,num_points,num_dim,num_dim);
603  auto cell_jac_det = af.buildStaticArray<Scalar,Cell,IP>("cell_jac_det",1,num_points);
604 
605  auto cell_basis_vector = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("cell_basis_vector",1,num_basis,num_points,num_dim);
606  auto cell_curl_basis_scalar = af.buildStaticArray<Scalar,Cell,BASIS,IP>("cell_curl_basis_scalar",1,num_basis,num_points);
607  auto cell_curl_basis_vector = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("cell_curl_basis_vector",1,num_basis,num_points,num_dim);
608 
609  auto cell_curl_basis_ref = af.buildArray<Scalar,BASIS,IP,Dim>("cell_curl_basis_ref",num_basis,num_points,num_dim);
610  auto cell_curl_basis_ref_scalar = af.buildStaticArray<Scalar,BASIS,IP>("cell_curl_basis_ref_scalar",num_basis,num_points);
611  auto cell_basis_ref_vector = af.buildArray<Scalar,BASIS,IP,Dim>("cell_basis_ref_vector",num_basis,num_points,num_dim);
612 
613  for(int cell=0;cell<num_cells;++cell){
614 
615  // =============================================
616  // Load external into cell-local arrays
617 
618  for(int p=0;p<num_points;++p)
619  for(int d=0;d<num_dim;++d)
620  for(int d2=0;d2<num_dim;++d2)
621  cell_jac(0,p,d,d2)=jac(cell,p,d,d2);
622  for(int p=0;p<num_points;++p)
623  for(int d=0;d<num_dim;++d)
624  for(int d2=0;d2<num_dim;++d2)
625  cell_jac_inv(0,p,d,d2)=jac_inv(cell,p,d,d2);
626  for(int p=0;p<num_points;++p)
627  cell_jac_det(0,p)=jac_det(cell,p);
628  for(int p=0;p<num_points;++p)
629  for(int d=0;d<num_dim;++d)
630  cell_cub_points(p,d)=cub_points(cell,p,d);
631 
632  // =============================================
633  // Load Reference Values
634 
635  intrepid_basis->getValues(cell_basis_ref_vector.get_view(),cell_cub_points.get_view(),Intrepid2::OPERATOR_VALUE);
636 
637  if(compute_derivatives){
638  if(num_dim==2){
639  intrepid_basis->getValues(cell_curl_basis_ref_scalar.get_view(),cell_cub_points.get_view(),Intrepid2::OPERATOR_CURL);
640  } else if(num_dim==3){
641  intrepid_basis->getValues(cell_curl_basis_ref.get_view(),cell_cub_points.get_view(),Intrepid2::OPERATOR_CURL);
642  }
643  }
644 
645  // =============================================
646  // Transform reference values to physical values
647 
648  fst::HCURLtransformVALUE(cell_basis_vector.get_view(),cell_jac_inv.get_view(),cell_basis_ref_vector.get_view());
649  for(int b=0;b<num_basis;++b)
650  for(int p=0;p<num_points;++p)
651  for(int d=0;d<num_dim;++d)
652  basis_vector(cell,b,p,d)=cell_basis_vector(0,b,p,d);
653 
654  if(compute_derivatives){
655  if(num_dim==2){
656  // note only volume deformation is needed!
657  // this relates directly to this being in
658  // the divergence space in 2D!
659  fst::HDIVtransformDIV(cell_curl_basis_scalar.get_view(),cell_jac_det.get_view(),cell_curl_basis_ref_scalar.get_view());
660  for(int b=0;b<num_basis;++b)
661  for(int p=0;p<num_points;++p)
662  curl_basis_scalar(cell,b,p)=cell_curl_basis_scalar(0,b,p);
663  } else if(num_dim==3) {
664  fst::HCURLtransformCURL(cell_curl_basis_vector.get_view(),cell_jac.get_view(),cell_jac_det.get_view(),cell_curl_basis_ref.get_view());
665  for(int b=0;b<num_basis;++b)
666  for(int p=0;p<num_points;++p)
667  for(int d=0;d<num_dim;++d)
668  curl_basis_vector(cell,b,p,d)=cell_curl_basis_vector(0,b,p,d);
669  } else {
670  TEUCHOS_TEST_FOR_EXCEPT_MSG(true,"panzer::BasisValues2::evaluateValues_HCurl : HCurl only setup for 2D and 3D.");
671  }
672  }
673  }
674 
675  if(build_weighted){
676  const std::pair<int,int> cell_range(0,num_cells);
677  auto s_weighted_basis_vector = Kokkos::subview(weighted_basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
678  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
679  auto s_basis_vector = Kokkos::subview(basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
680  fst::multiplyMeasure(s_weighted_basis_vector,s_weighted_measure,s_basis_vector);
681  if(compute_derivatives){
682  if(num_dim==2){
683  auto s_weighted_curl_basis_scalar = Kokkos::subview(weighted_curl_basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
684  auto s_curl_basis_scalar = Kokkos::subview(curl_basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
685  fst::multiplyMeasure(s_weighted_curl_basis_scalar,s_weighted_measure,s_curl_basis_scalar);
686  } else if(num_dim==3){
687  auto s_weighted_curl_basis_vector = Kokkos::subview(weighted_curl_basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
688  auto s_curl_basis_vector = Kokkos::subview(curl_basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
689  fst::multiplyMeasure(s_weighted_curl_basis_vector,s_weighted_measure,s_curl_basis_vector);
690  }
691  }
692  }
693 
694 }
695 
696 template <typename Scalar>
698 evaluateValues_HDiv(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
699  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
700  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
701  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
702  const int in_num_cells)
703 {
704 
705  TEUCHOS_ASSERT(getElementSpace() == PureBasis::HDIV);
706 
707  typedef Intrepid2::FunctionSpaceTools<PHX::Device::execution_space> fst;
708  MDFieldArrayFactory af("",ddims_,true);
709 
710  const panzer::PureBasis & basis = *(basis_layout->getBasis());
711 
712  const int num_points = basis_layout->numPoints();
713  const int num_basis = basis.cardinality();
714  const int num_dim = basis_layout->dimension();
715  const int num_cells = in_num_cells < 0 ? basis_layout->numCells() : in_num_cells;
716 
717  auto cell_cub_points = af.buildStaticArray<Scalar,IP,Dim>("cell_cub_points",num_points,num_dim);
718  auto cell_jac = af.buildStaticArray<Scalar,Cell,IP,Dim,Dim>("cell_jac",1,num_points,num_dim,num_dim);
719  auto cell_jac_det = af.buildStaticArray<Scalar,Cell,IP>("cell_jac_det",1,num_points);
720 
721  auto cell_basis_vector = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("cell_basis_vector",1,num_basis,num_points,num_dim);
722  auto cell_div_basis = af.buildStaticArray<Scalar,Cell,BASIS,IP>("cell_div_basis",1,num_basis,num_points);
723 
724  auto cell_basis_ref_vector = af.buildArray<Scalar,BASIS,IP,Dim>("cell_basis_ref_vector",num_basis,num_points,num_dim);
725  auto cell_div_basis_ref = af.buildStaticArray<Scalar,BASIS,IP>("cell_div_basis_ref",num_basis,num_points);
726 
727  for(int cell=0;cell<num_cells;++cell){
728 
729  // =============================================
730  // Load external into cell-local arrays
731 
732  for(int p=0;p<num_points;++p)
733  for(int d=0;d<num_dim;++d)
734  for(int d2=0;d2<num_dim;++d2)
735  cell_jac(0,p,d,d2)=jac(cell,p,d,d2);
736  for(int p=0;p<num_points;++p)
737  cell_jac_det(0,p)=jac_det(cell,p);
738  for(int p=0;p<num_points;++p)
739  for(int d=0;d<num_dim;++d)
740  cell_cub_points(p,d)=cub_points(cell,p,d);
741  // =============================================
742  // Load Reference Values
743 
744  intrepid_basis->getValues(cell_basis_ref_vector.get_view(),cell_cub_points.get_view(),Intrepid2::OPERATOR_VALUE);
745 
746  if(compute_derivatives){
747  intrepid_basis->getValues(cell_div_basis_ref.get_view(),cell_cub_points.get_view(),Intrepid2::OPERATOR_DIV);
748  }
749 
750  // =============================================
751  // Transform reference values to physical values
752 
753  fst::HDIVtransformVALUE(cell_basis_vector.get_view(),cell_jac.get_view(),cell_jac_det.get_view(),cell_basis_ref_vector.get_view());
754  for(int b=0;b<num_basis;++b)
755  for(int p=0;p<num_points;++p)
756  for(int d=0;d<num_dim;++d)
757  basis_vector(cell,b,p,d)=cell_basis_vector(0,b,p,d);
758 
759  if(compute_derivatives){
760  fst::HDIVtransformDIV(cell_div_basis.get_view(),cell_jac_det.get_view(),cell_div_basis_ref.get_view());
761  for(int b=0;b<num_basis;++b)
762  for(int p=0;p<num_points;++p)
763  div_basis(cell,b,p)=cell_div_basis(0,b,p);
764  }
765  }
766 
767  if(build_weighted){
768  const std::pair<int,int> cell_range(0,num_cells);
769  auto s_weighted_basis_vector = Kokkos::subview(weighted_basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
770  auto s_weighted_measure = Kokkos::subview(weighted_measure.get_view(),cell_range,Kokkos::ALL());
771  auto s_basis_vector = Kokkos::subview(basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
772  fst::multiplyMeasure(s_weighted_basis_vector,s_weighted_measure,s_basis_vector);
773  if(compute_derivatives){
774  auto s_weighted_div_basis = Kokkos::subview(weighted_div_basis.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
775  auto s_div_basis = Kokkos::subview(div_basis.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
776  fst::multiplyMeasure(s_weighted_div_basis,s_weighted_measure,s_div_basis);
777  }
778  }
779 
780 }
781 
782 
783 
784 
785 
786 
787 
788 
789 
790 
791 
792 
793 
794 
795 
796 
797 
798 
799 
800 
801 
802 
803 
804 
805 
806 
807 
808 
809 
810 
811 
812 
813 
814 
815 
816 
817 
818 
819 
820 
821 
822 
823 
824 template <typename Scalar>
826 evaluateValuesCV(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cell_cub_points,
827  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
828  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
829  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv)
830 {
831  MDFieldArrayFactory af("",ddims_,true);
832 
833  int num_ip = basis_layout->numPoints();
834  int num_card = basis_layout->cardinality();
835  int num_dim = basis_layout->dimension();
836 
837  size_type num_cells = jac.extent(0);
838 
839  PureBasis::EElementSpace elmtspace = getElementSpace();
840  ArrayDynamic dyn_cub_points = af.buildArray<Scalar,IP,Dim>("dyn_cub_points", num_ip, num_dim);
841 
842  // Integration points are located on physical cells rather than reference cells,
843  // so we evaluate the basis in a loop over cells.
844  for (size_type icell = 0; icell < num_cells; ++icell)
845  {
846  for (int ip = 0; ip < num_ip; ++ip)
847  for (int d = 0; d < num_dim; ++d)
848  dyn_cub_points(ip,d) = cell_cub_points(icell,ip,d);
849 
850  if(elmtspace==PureBasis::CONST) {
851  ArrayDynamic dyn_basis_ref_scalar = af.buildArray<Scalar,BASIS,IP>("dyn_basis_ref_scalar",num_card,num_ip);
852 
853  intrepid_basis->getValues(dyn_basis_ref_scalar.get_view(),
854  dyn_cub_points.get_view(),
855  Intrepid2::OPERATOR_VALUE);
856 
857  // transform values method just transfers values to array with cell index - no need to call
858  for (int b = 0; b < num_card; ++b)
859  for (int ip = 0; ip < num_ip; ++ip)
860  basis_scalar(icell,b,ip) = dyn_basis_ref_scalar(b,ip);
861 
862  }
863  if(elmtspace==PureBasis::HVOL) {
864  ArrayDynamic dyn_basis_ref_scalar = af.buildArray<Scalar,BASIS,IP>("dyn_basis_ref_scalar",num_card,num_ip);
865 
866  intrepid_basis->getValues(dyn_basis_ref_scalar.get_view(),
867  dyn_cub_points.get_view(),
868  Intrepid2::OPERATOR_VALUE);
869 
870  int one_cell= 1;
871  ArrayDynamic dyn_basis_scalar = af.buildArray<Scalar,Cell,BASIS,IP>("dyn_basis_vector",one_cell,num_card,num_ip);
872  ArrayDynamic dyn_jac = af.buildArray<Scalar,Cell,IP,Dim,Dim>("dyn_jac",one_cell,num_ip,num_dim,num_dim);
873  ArrayDynamic dyn_jac_det = af.buildArray<Scalar,Cell,IP>("dyn_jac_det",one_cell,num_ip);
874 
875  int cellInd = 0;
876  for (int ip = 0; ip < num_ip; ++ip)
877  dyn_jac_det(cellInd,ip) = jac_det(icell,ip);
878 
879  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::HVOLtransformVALUE(dyn_basis_scalar.get_view(),
880  dyn_jac_det.get_view(),
881  dyn_basis_ref_scalar.get_view());
882 
883  for (int b = 0; b < num_card; ++b)
884  for (int ip = 0; ip < num_ip; ++ip)
885  basis_scalar(icell,b,ip) = dyn_basis_scalar(0,b,ip);
886 
887  }
888  if(elmtspace==PureBasis::HGRAD) {
889  ArrayDynamic dyn_basis_ref_scalar = af.buildArray<Scalar,BASIS,IP>("dyn_basis_ref_scalar",num_card,num_ip);
890 
891  intrepid_basis->getValues(dyn_basis_ref_scalar.get_view(),
892  dyn_cub_points.get_view(),
893  Intrepid2::OPERATOR_VALUE);
894 
895  // transform values method just transfers values to array with cell index - no need to call
896  for (int b = 0; b < num_card; ++b)
897  for (int ip = 0; ip < num_ip; ++ip)
898  basis_scalar(icell,b,ip) = dyn_basis_ref_scalar(b,ip);
899 
900  if(compute_derivatives) {
901 
902  int one_cell = 1;
903  ArrayDynamic dyn_grad_basis_ref = af.buildArray<Scalar,BASIS,IP,Dim>("dyn_grad_basis_ref",num_card,num_ip,num_dim);
904  ArrayDynamic dyn_grad_basis = af.buildArray<Scalar,Cell,BASIS,IP,Dim>("dyn_grad_basis",one_cell,num_card,num_ip,num_dim);
905  ArrayDynamic dyn_jac_inv = af.buildArray<Scalar,Cell,IP,Dim,Dim>("dyn_jac_inv",one_cell,num_ip,num_dim,num_dim);
906 
907  intrepid_basis->getValues(dyn_grad_basis_ref.get_view(),
908  dyn_cub_points.get_view(),
909  Intrepid2::OPERATOR_GRAD);
910 
911  int cellInd = 0;
912  for (int ip = 0; ip < num_ip; ++ip)
913  for (int d1 = 0; d1 < num_dim; ++d1)
914  for (int d2 = 0; d2 < num_dim; ++d2)
915  dyn_jac_inv(cellInd,ip,d1,d2) = jac_inv(icell,ip,d1,d2);
916 
917  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::HGRADtransformGRAD<Scalar>(dyn_grad_basis.get_view(),
918  dyn_jac_inv.get_view(),
919  dyn_grad_basis_ref.get_view());
920 
921  for (int b = 0; b < num_card; ++b)
922  for (int ip = 0; ip < num_ip; ++ip)
923  for (int d = 0; d < num_dim; ++d)
924  grad_basis(icell,b,ip,d) = dyn_grad_basis(0,b,ip,d);
925 
926  }
927  }
928  else if(elmtspace==PureBasis::HCURL) {
929  ArrayDynamic dyn_basis_ref_vector = af.buildArray<Scalar,BASIS,IP,Dim>("dyn_basis_ref_vector",num_card,num_ip,num_dim);
930 
931  intrepid_basis->getValues(dyn_basis_ref_vector.get_view(),
932  dyn_cub_points.get_view(),
933  Intrepid2::OPERATOR_VALUE);
934 
935  const int one_cell = 1;
936  ArrayDynamic dyn_basis_vector = af.buildArray<Scalar,Cell,BASIS,IP,Dim>("dyn_basis_vector",one_cell,num_card,num_ip,num_dim);
937  ArrayDynamic dyn_jac_inv = af.buildArray<Scalar,Cell,IP,Dim,Dim>("dyn_jac_inv",one_cell,num_ip,num_dim,num_dim);
938 
939  const int cellInd = 0;
940  for (int ip = 0; ip < num_ip; ++ip)
941  for (int d1 = 0; d1 < num_dim; ++d1)
942  for (int d2 = 0; d2 < num_dim; ++d2)
943  dyn_jac_inv(cellInd,ip,d1,d2) = jac_inv(icell,ip,d1,d2);
944 
945  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::HCURLtransformVALUE(dyn_basis_vector.get_view(),
946  dyn_jac_inv.get_view(),
947  dyn_basis_ref_vector.get_view());
948 
949  for (int b = 0; b < num_card; ++b)
950  for (int ip = 0; ip < num_ip; ++ip)
951  for (int d = 0; d < num_dim; ++d)
952  basis_vector(icell,b,ip,d) = dyn_basis_vector(0,b,ip,d);
953 
954  if(compute_derivatives && num_dim ==2) {
955 
956  ArrayDynamic dyn_curl_basis_ref_scalar = af.buildArray<Scalar,BASIS,IP>("dyn_curl_basis_ref_scalar",num_card,num_ip);
957  ArrayDynamic dyn_curl_basis_scalar = af.buildArray<Scalar,Cell,BASIS,IP>("dyn_curl_basis_scalar",one_cell,num_card,num_ip);
958  ArrayDynamic dyn_jac_det = af.buildArray<Scalar,Cell,IP>("dyn_jac_det",one_cell,num_ip);
959 
960  intrepid_basis->getValues(dyn_curl_basis_ref_scalar.get_view(),
961  dyn_cub_points.get_view(),
962  Intrepid2::OPERATOR_CURL);
963 
964  for (int ip = 0; ip < num_ip; ++ip)
965  dyn_jac_det(cellInd,ip) = jac_det(icell,ip);
966 
967  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::HDIVtransformDIV(dyn_curl_basis_scalar.get_view(),
968  dyn_jac_det.get_view(),
969  dyn_curl_basis_ref_scalar.get_view());
970 
971  for (int b = 0; b < num_card; ++b)
972  for (int ip = 0; ip < num_ip; ++ip)
973  curl_basis_scalar(icell,b,ip) = dyn_curl_basis_scalar(0,b,ip);
974 
975  }
976  if(compute_derivatives && num_dim ==3) {
977 
978  ArrayDynamic dyn_curl_basis_ref = af.buildArray<Scalar,BASIS,IP,Dim>("dyn_curl_basis_ref_vector",num_card,num_ip,num_dim);
979  ArrayDynamic dyn_curl_basis = af.buildArray<Scalar,Cell,BASIS,IP,Dim>("dyn_curl_basis_vector",one_cell,num_card,num_ip,num_dim);
980  ArrayDynamic dyn_jac_det = af.buildArray<Scalar,Cell,IP>("dyn_jac_det",one_cell,num_ip);
981  ArrayDynamic dyn_jac = af.buildArray<Scalar,Cell,IP,Dim,Dim>("dyn_jac",one_cell,num_ip,num_dim,num_dim);
982 
983  intrepid_basis->getValues(dyn_curl_basis_ref.get_view(),
984  dyn_cub_points.get_view(),
985  Intrepid2::OPERATOR_CURL);
986 
987  for (int ip = 0; ip < num_ip; ++ip)
988  {
989  dyn_jac_det(cellInd,ip) = jac_det(icell,ip);
990  for (int d1 = 0; d1 < num_dim; ++d1)
991  for (int d2 = 0; d2 < num_dim; ++d2)
992  dyn_jac(cellInd,ip,d1,d2) = jac(icell,ip,d1,d2);
993  }
994 
995  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::HCURLtransformCURL(dyn_curl_basis.get_view(),
996  dyn_jac.get_view(),
997  dyn_jac_det.get_view(),
998  dyn_curl_basis_ref.get_view());
999 
1000  for (int b = 0; b < num_card; ++b)
1001  for (int ip = 0; ip < num_ip; ++ip)
1002  for (int d = 0; d < num_dim; ++d)
1003  curl_basis_vector(icell,b,ip,d) = dyn_curl_basis(0,b,ip,d);
1004 
1005  }
1006 
1007  }
1008  else if(elmtspace==PureBasis::HDIV) {
1009 
1010  ArrayDynamic dyn_basis_ref_vector = af.buildArray<Scalar,BASIS,IP,Dim>("dyn_basis_ref_vector",num_card,num_ip,num_dim);
1011 
1012  intrepid_basis->getValues(dyn_basis_ref_vector.get_view(),
1013  dyn_cub_points.get_view(),
1014  Intrepid2::OPERATOR_VALUE);
1015 
1016  int one_cell= 1;
1017  ArrayDynamic dyn_basis_vector = af.buildArray<Scalar,Cell,BASIS,IP,Dim>("dyn_basis_vector",one_cell,num_card,num_ip,num_dim);
1018  ArrayDynamic dyn_jac = af.buildArray<Scalar,Cell,IP,Dim,Dim>("dyn_jac",one_cell,num_ip,num_dim,num_dim);
1019  ArrayDynamic dyn_jac_det = af.buildArray<Scalar,Cell,IP>("dyn_jac_det",one_cell,num_ip);
1020 
1021  int cellInd = 0;
1022  for (int ip = 0; ip < num_ip; ++ip)
1023  {
1024  dyn_jac_det(cellInd,ip) = jac_det(icell,ip);
1025  for (int d1 = 0; d1 < num_dim; ++d1)
1026  for (int d2 = 0; d2 < num_dim; ++d2)
1027  dyn_jac(cellInd,ip,d1,d2) = jac(icell,ip,d1,d2);
1028  }
1029 
1030  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::HDIVtransformVALUE(dyn_basis_vector.get_view(),
1031  dyn_jac.get_view(),
1032  dyn_jac_det.get_view(),
1033  dyn_basis_ref_vector.get_view());
1034 
1035  for (int b = 0; b < num_card; ++b)
1036  for (int ip = 0; ip < num_ip; ++ip)
1037  for (int d = 0; d < num_dim; ++d)
1038  basis_vector(icell,b,ip,d) = dyn_basis_vector(0,b,ip,d);
1039 
1040  if(compute_derivatives) {
1041 
1042  ArrayDynamic dyn_div_basis_ref = af.buildArray<Scalar,BASIS,IP>("dyn_div_basis_ref_scalar",num_card,num_ip);
1043  ArrayDynamic dyn_div_basis = af.buildArray<Scalar,Cell,BASIS,IP>("dyn_div_basis_scalar",one_cell,num_card,num_ip);
1044 
1045  intrepid_basis->getValues(dyn_div_basis_ref.get_view(),
1046  dyn_cub_points.get_view(),
1047  Intrepid2::OPERATOR_DIV);
1048 
1049  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::HDIVtransformDIV<Scalar>(dyn_div_basis.get_view(),
1050  dyn_jac_det.get_view(),
1051  dyn_div_basis_ref.get_view());
1052 
1053  for (int b = 0; b < num_card; ++b)
1054  for (int ip = 0; ip < num_ip; ++ip)
1055  div_basis(icell,b,ip) = dyn_div_basis(0,b,ip);
1056 
1057  }
1058 
1059  }
1060  else { TEUCHOS_ASSERT(false); }
1061 
1062  } // cell loop
1063 
1064 }
1065 
1066 template <typename Scalar>
1068 evaluateReferenceValues(const PHX::MDField<Scalar,IP,Dim> & cub_points,bool in_compute_derivatives,bool use_vertex_coordinates)
1069 {
1070  MDFieldArrayFactory af("",ddims_,true);
1071 
1072  int num_quad = basis_layout->numPoints();
1073  int num_dim = basis_layout->dimension();
1074  int num_card = basis_layout->cardinality();
1075 
1076  ArrayDynamic dyn_cub_points = af.buildArray<Scalar,IP,Dim>("dyn_cub_points", num_quad,num_dim);
1077 
1078  for (int ip = 0; ip < num_quad; ++ip)
1079  for (int d = 0; d < num_dim; ++d)
1080  dyn_cub_points(ip,d) = cub_points(ip,d);
1081 
1082  PureBasis::EElementSpace elmtspace = getElementSpace();
1083  if(elmtspace==PureBasis::HGRAD || elmtspace==PureBasis::CONST || elmtspace==PureBasis::HVOL) {
1084  ArrayDynamic dyn_basis_ref_scalar = af.buildArray<Scalar,BASIS,IP>("dyn_basis_ref_scalar",num_card,num_quad);
1085 
1086  intrepid_basis->getValues(dyn_basis_ref_scalar.get_view(),
1087  dyn_cub_points.get_view(),
1088  Intrepid2::OPERATOR_VALUE);
1089 
1090  for (int b = 0; b < num_card; ++b)
1091  for (int ip = 0; ip < num_quad; ++ip)
1092  basis_ref_scalar(b,ip) = dyn_basis_ref_scalar(b,ip);
1093  }
1094  else if(elmtspace==PureBasis::HDIV || elmtspace==PureBasis::HCURL) {
1095  ArrayDynamic dyn_basis_ref_vector = af.buildArray<Scalar,BASIS,IP,Dim>("dyn_basis_ref_vector",num_card,num_quad,num_dim);
1096 
1097  intrepid_basis->getValues(dyn_basis_ref_vector.get_view(),
1098  dyn_cub_points.get_view(),
1099  Intrepid2::OPERATOR_VALUE);
1100 
1101  for (int b = 0; b < num_card; ++b)
1102  for (int ip = 0; ip < num_quad; ++ip)
1103  for (int d = 0; d < num_dim; ++d)
1104  basis_ref_vector(b,ip,d) = dyn_basis_ref_vector(b,ip,d);
1105  }
1106  else { TEUCHOS_ASSERT(false); }
1107 
1108  if(elmtspace==PureBasis::HGRAD && in_compute_derivatives) {
1109  ArrayDynamic dyn_grad_basis_ref = af.buildArray<Scalar,BASIS,IP,Dim>("dyn_basis_ref_vector",num_card,num_quad,num_dim);
1110 
1111  intrepid_basis->getValues(dyn_grad_basis_ref.get_view(),
1112  dyn_cub_points.get_view(),
1113  Intrepid2::OPERATOR_GRAD);
1114 
1115  for (int b = 0; b < num_card; ++b)
1116  for (int ip = 0; ip < num_quad; ++ip)
1117  for (int d = 0; d < num_dim; ++d)
1118  grad_basis_ref(b,ip,d) = dyn_grad_basis_ref(b,ip,d);
1119  }
1120  else if(elmtspace==PureBasis::HCURL && in_compute_derivatives && num_dim==2) {
1121  ArrayDynamic dyn_curl_basis_ref = af.buildArray<Scalar,BASIS,IP>("dyn_curl_basis_ref_scalar",num_card,num_quad);
1122 
1123  intrepid_basis->getValues(dyn_curl_basis_ref.get_view(),
1124  dyn_cub_points.get_view(),
1125  Intrepid2::OPERATOR_CURL);
1126 
1127  for (int b = 0; b < num_card; ++b)
1128  for (int ip = 0; ip < num_quad; ++ip)
1129  curl_basis_ref_scalar(b,ip) = dyn_curl_basis_ref(b,ip);
1130  }
1131  else if(elmtspace==PureBasis::HCURL && in_compute_derivatives && num_dim==3) {
1132  ArrayDynamic dyn_curl_basis_ref = af.buildArray<Scalar,BASIS,IP,Dim>("dyn_curl_basis_ref_vector",num_card,num_quad,num_dim);
1133 
1134  intrepid_basis->getValues(dyn_curl_basis_ref.get_view(),
1135  dyn_cub_points.get_view(),
1136  Intrepid2::OPERATOR_CURL);
1137 
1138  for (int b = 0; b < num_card; ++b)
1139  for (int ip = 0; ip < num_quad; ++ip)
1140  for (int d = 0; d < num_dim; ++d)
1141  curl_basis_ref_vector(b,ip,d) = dyn_curl_basis_ref(b,ip,d);
1142  }
1143  else if(elmtspace==PureBasis::HDIV && in_compute_derivatives) {
1144  ArrayDynamic dyn_div_basis_ref = af.buildArray<Scalar,BASIS,IP>("dyn_div_basis_ref_scalar",num_card,num_quad);
1145 
1146  intrepid_basis->getValues(dyn_div_basis_ref.get_view(),
1147  dyn_cub_points.get_view(),
1148  Intrepid2::OPERATOR_DIV);
1149 
1150  for (int b = 0; b < num_card; ++b)
1151  for (int ip = 0; ip < num_quad; ++ip)
1152  div_basis_ref(b,ip) = dyn_div_basis_ref(b,ip);
1153  }
1154 
1155 
1156  if(use_vertex_coordinates) {
1157  // Intrepid removes fad types from the coordinate scalar type. We
1158  // pull the actual field scalar type from the basis object to be
1159  // consistent.
1160  if (elmtspace != PureBasis::CONST) {
1161  using coordsScalarType = typename Intrepid2::Basis<PHX::Device::execution_space,Scalar,Scalar>::scalarType;
1162  auto dyn_basis_coordinates_ref = af.buildArray<coordsScalarType,BASIS,Dim>("basis_coordinates_ref",
1163  basis_coordinates_ref.extent(0),
1164  basis_coordinates_ref.extent(1));
1165  intrepid_basis->getDofCoords(dyn_basis_coordinates_ref.get_view());
1166 
1167  // fill in basis coordinates
1168  for (int i = 0; i < basis_coordinates_ref.extent_int(0); ++i)
1169  for (int j = 0; j < basis_coordinates_ref.extent_int(1); ++j)
1170  basis_coordinates_ref(i,j) = dyn_basis_coordinates_ref(i,j);
1171  }
1172  }
1173 
1174  references_evaluated = true;
1175 }
1176 
1177 // method for applying orientations
1178 template <typename Scalar>
1180 applyOrientations(const std::vector<Intrepid2::Orientation> & orientations,
1181  const int in_num_cells)
1182 {
1183  if (!intrepid_basis->requireOrientation())
1184  return;
1185 
1186  typedef Intrepid2::OrientationTools<PHX::Device> ots;
1187  const PureBasis::EElementSpace elmtspace = getElementSpace();
1188 
1189  // orientation (right now std vector) is created using push_back method.
1190  // thus, its size is the actual number of elements to be applied.
1191  // on the other hand, basis_layout num cell indicates workset size.
1192  // to get right size of cells, use minimum of them.
1193  const int num_cell_basis_layout = in_num_cells < 0 ? basis_layout->numCells() : in_num_cells;
1194  const int num_cell_orientation = orientations.size();
1195  const int num_cell = num_cell_basis_layout < num_cell_orientation ? num_cell_basis_layout : num_cell_orientation;
1196  const int num_dim = basis_layout->dimension();
1197  const Kokkos::pair<int,int> range_cell(0, num_cell);
1198 
1199  // Kokkos::DynRankView<Intrepid2::Orientation,PHX::Device>
1200  // drv_orts((Intrepid2::Orientation*)orientations.data(), num_cell);
1201  Kokkos::DynRankView<Intrepid2::Orientation,PHX::Device> drv_orts("drv_orts", num_cell);
1202  auto host_drv_orts = Kokkos::create_mirror_view(drv_orts);
1203  for (size_t i=0; i < drv_orts.size(); ++i)
1204  host_drv_orts(i) = orientations[i];
1205  Kokkos::deep_copy(drv_orts,host_drv_orts);
1206  PHX::Device::fence();
1207 
1211  if (elmtspace==PureBasis::HGRAD) {
1212  {
1213  {
1214  auto drv_basis_scalar = Kokkos::subview(basis_scalar.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL());
1215  auto drv_basis_scalar_tmp = Kokkos::createDynRankView(basis_scalar.get_view(),
1216  "drv_basis_scalar_tmp",
1217  drv_basis_scalar.extent(0), // C
1218  drv_basis_scalar.extent(1), // F
1219  drv_basis_scalar.extent(2)); // P
1220  Kokkos::deep_copy(drv_basis_scalar_tmp, drv_basis_scalar);
1221  ots::modifyBasisByOrientation(drv_basis_scalar,
1222  drv_basis_scalar_tmp,
1223  drv_orts,
1224  intrepid_basis);
1225  }
1226  if(build_weighted) {
1227  auto drv_basis_scalar = Kokkos::subview(weighted_basis_scalar.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL());
1228  auto drv_basis_scalar_tmp = Kokkos::createDynRankView(weighted_basis_scalar.get_view(),
1229  "drv_basis_scalar_tmp",
1230  drv_basis_scalar.extent(0), // C
1231  drv_basis_scalar.extent(1), // F
1232  drv_basis_scalar.extent(2)); // P
1233  Kokkos::deep_copy(drv_basis_scalar_tmp, drv_basis_scalar);
1234  ots::modifyBasisByOrientation(drv_basis_scalar,
1235  drv_basis_scalar_tmp,
1236  drv_orts,
1237  intrepid_basis);
1238  }
1239 
1240  }
1241 
1242  if (compute_derivatives) {
1243  {
1244  auto drv_grad_basis = Kokkos::subview(grad_basis.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1245  auto drv_grad_basis_tmp = Kokkos::createDynRankView(grad_basis.get_view(),
1246  "drv_grad_basis_tmp",
1247  drv_grad_basis.extent(0), // C
1248  drv_grad_basis.extent(1), // F
1249  drv_grad_basis.extent(2), // P
1250  drv_grad_basis.extent(3)); // D
1251  Kokkos::deep_copy(drv_grad_basis_tmp, drv_grad_basis);
1252  ots::modifyBasisByOrientation(drv_grad_basis,
1253  drv_grad_basis_tmp,
1254  drv_orts,
1255  intrepid_basis);
1256  }
1257  if(build_weighted) {
1258  auto drv_grad_basis = Kokkos::subview(weighted_grad_basis.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1259  auto drv_grad_basis_tmp = Kokkos::createDynRankView(weighted_grad_basis.get_view(),
1260  "drv_grad_basis_tmp",
1261  drv_grad_basis.extent(0), // C
1262  drv_grad_basis.extent(1), // F
1263  drv_grad_basis.extent(2), // P
1264  drv_grad_basis.extent(3)); // D
1265  Kokkos::deep_copy(drv_grad_basis_tmp, drv_grad_basis);
1266  ots::modifyBasisByOrientation(drv_grad_basis,
1267  drv_grad_basis_tmp,
1268  drv_orts,
1269  intrepid_basis);
1270  }
1271  }
1272  }
1273 
1277  else if (elmtspace==PureBasis::HCURL && num_dim==2) {
1278  {
1279  {
1280  auto drv_basis_vector = Kokkos::subview(basis_vector.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1281  auto drv_basis_vector_tmp = Kokkos::createDynRankView(basis_vector.get_view(),
1282  "drv_basis_vector_tmp",
1283  drv_basis_vector.extent(0), // C
1284  drv_basis_vector.extent(1), // F
1285  drv_basis_vector.extent(2), // P
1286  drv_basis_vector.extent(3)); // D
1287  Kokkos::deep_copy(drv_basis_vector_tmp, drv_basis_vector);
1288  ots::modifyBasisByOrientation(drv_basis_vector,
1289  drv_basis_vector_tmp,
1290  drv_orts,
1291  intrepid_basis);
1292  }
1293  if(build_weighted) {
1294  auto drv_basis_vector = Kokkos::subview(weighted_basis_vector.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1295  auto drv_basis_vector_tmp = Kokkos::createDynRankView(weighted_basis_vector.get_view(),
1296  "drv_basis_vector_tmp",
1297  drv_basis_vector.extent(0), // C
1298  drv_basis_vector.extent(1), // F
1299  drv_basis_vector.extent(2), // P
1300  drv_basis_vector.extent(3)); // D
1301  Kokkos::deep_copy(drv_basis_vector_tmp, drv_basis_vector);
1302  ots::modifyBasisByOrientation(drv_basis_vector,
1303  drv_basis_vector_tmp,
1304  drv_orts,
1305  intrepid_basis);
1306  }
1307  }
1308 
1309  if (compute_derivatives) {
1310  {
1311  auto drv_curl_basis_scalar = Kokkos::subview(curl_basis_scalar.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL());
1312  auto drv_curl_basis_scalar_tmp = Kokkos::createDynRankView(curl_basis_scalar.get_view(),
1313  "drv_curl_basis_scalar_tmp",
1314  drv_curl_basis_scalar.extent(0), // C
1315  drv_curl_basis_scalar.extent(1), // F
1316  drv_curl_basis_scalar.extent(2)); // P
1317  Kokkos::deep_copy(drv_curl_basis_scalar_tmp, drv_curl_basis_scalar);
1318  ots::modifyBasisByOrientation(drv_curl_basis_scalar,
1319  drv_curl_basis_scalar_tmp,
1320  drv_orts,
1321  intrepid_basis);
1322  }
1323 
1324  if(build_weighted) {
1325  auto drv_curl_basis_scalar = Kokkos::subview(weighted_curl_basis_scalar.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL());
1326  auto drv_curl_basis_scalar_tmp = Kokkos::createDynRankView(weighted_curl_basis_scalar.get_view(),
1327  "drv_curl_basis_scalar_tmp",
1328  drv_curl_basis_scalar.extent(0), // C
1329  drv_curl_basis_scalar.extent(1), // F
1330  drv_curl_basis_scalar.extent(2)); // P
1331  Kokkos::deep_copy(drv_curl_basis_scalar_tmp, drv_curl_basis_scalar);
1332  ots::modifyBasisByOrientation(drv_curl_basis_scalar,
1333  drv_curl_basis_scalar_tmp,
1334  drv_orts,
1335  intrepid_basis);
1336  }
1337  }
1338  }
1339 
1343  else if (elmtspace==PureBasis::HCURL && num_dim==3) {
1344  {
1345  {
1346  auto drv_basis_vector = Kokkos::subview(basis_vector.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1347  auto drv_basis_vector_tmp = Kokkos::createDynRankView(basis_vector.get_view(),
1348  "drv_basis_vector_tmp",
1349  drv_basis_vector.extent(0), // C
1350  drv_basis_vector.extent(1), // F
1351  drv_basis_vector.extent(2), // P
1352  drv_basis_vector.extent(3)); // D
1353  Kokkos::deep_copy(drv_basis_vector_tmp, drv_basis_vector);
1354  ots::modifyBasisByOrientation(drv_basis_vector,
1355  drv_basis_vector_tmp,
1356  drv_orts,
1357  intrepid_basis);
1358  }
1359  if(build_weighted) {
1360  auto drv_basis_vector = Kokkos::subview(weighted_basis_vector.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1361  auto drv_basis_vector_tmp = Kokkos::createDynRankView(weighted_basis_vector.get_view(),
1362  "drv_basis_vector_tmp",
1363  drv_basis_vector.extent(0), // C
1364  drv_basis_vector.extent(1), // F
1365  drv_basis_vector.extent(2), // P
1366  drv_basis_vector.extent(3)); // D
1367  Kokkos::deep_copy(drv_basis_vector_tmp, drv_basis_vector);
1368  ots::modifyBasisByOrientation(drv_basis_vector,
1369  drv_basis_vector_tmp,
1370  drv_orts,
1371  intrepid_basis);
1372  }
1373  }
1374 
1375  if (compute_derivatives) {
1376  {
1377  auto drv_curl_basis_vector = Kokkos::subview(curl_basis_vector.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1378  auto drv_curl_basis_vector_tmp = Kokkos::createDynRankView(curl_basis_vector.get_view(),
1379  "drv_curl_basis_vector_tmp",
1380  drv_curl_basis_vector.extent(0), // C
1381  drv_curl_basis_vector.extent(1), // F
1382  drv_curl_basis_vector.extent(2), // P
1383  drv_curl_basis_vector.extent(3)); // D
1384  Kokkos::deep_copy(drv_curl_basis_vector_tmp, drv_curl_basis_vector);
1385  ots::modifyBasisByOrientation(drv_curl_basis_vector,
1386  drv_curl_basis_vector_tmp,
1387  drv_orts,
1388  intrepid_basis);
1389  }
1390  if(build_weighted) {
1391  auto drv_curl_basis_vector = Kokkos::subview(weighted_curl_basis_vector.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1392  auto drv_curl_basis_vector_tmp = Kokkos::createDynRankView(weighted_curl_basis_vector.get_view(),
1393  "drv_curl_basis_vector_tmp",
1394  drv_curl_basis_vector.extent(0), // C
1395  drv_curl_basis_vector.extent(1), // F
1396  drv_curl_basis_vector.extent(2), // P
1397  drv_curl_basis_vector.extent(3)); // D
1398  Kokkos::deep_copy(drv_curl_basis_vector_tmp, drv_curl_basis_vector);
1399  ots::modifyBasisByOrientation(drv_curl_basis_vector,
1400  drv_curl_basis_vector_tmp,
1401  drv_orts,
1402  intrepid_basis);
1403  }
1404  }
1405  }
1409  else if (elmtspace==PureBasis::HDIV) {
1410  {
1411  {
1412  auto drv_basis_vector = Kokkos::subview(basis_vector.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1413  auto drv_basis_vector_tmp = Kokkos::createDynRankView(basis_vector.get_view(),
1414  "drv_basis_vector_tmp",
1415  drv_basis_vector.extent(0), // C
1416  drv_basis_vector.extent(1), // F
1417  drv_basis_vector.extent(2), // P
1418  drv_basis_vector.extent(3)); // D
1419  Kokkos::deep_copy(drv_basis_vector_tmp, drv_basis_vector);
1420  ots::modifyBasisByOrientation(drv_basis_vector,
1421  drv_basis_vector_tmp,
1422  drv_orts,
1423  intrepid_basis);
1424  }
1425  if(build_weighted) {
1426  auto drv_basis_vector = Kokkos::subview(weighted_basis_vector.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1427  auto drv_basis_vector_tmp = Kokkos::createDynRankView(weighted_basis_vector.get_view(),
1428  "drv_basis_vector_tmp",
1429  drv_basis_vector.extent(0), // C
1430  drv_basis_vector.extent(1), // F
1431  drv_basis_vector.extent(2), // P
1432  drv_basis_vector.extent(3)); // D
1433  Kokkos::deep_copy(drv_basis_vector_tmp, drv_basis_vector);
1434  ots::modifyBasisByOrientation(drv_basis_vector,
1435  drv_basis_vector_tmp,
1436  drv_orts,
1437  intrepid_basis);
1438  }
1439  }
1440  if (compute_derivatives) {
1441  {
1442  auto drv_div_basis = Kokkos::subview(div_basis.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL());
1443  auto drv_div_basis_tmp = Kokkos::createDynRankView(div_basis.get_view(),
1444  "drv_div_basis_tmp",
1445  drv_div_basis.extent(0), // C
1446  drv_div_basis.extent(1), // F
1447  drv_div_basis.extent(2)); // P
1448  Kokkos::deep_copy(drv_div_basis_tmp, drv_div_basis);
1449  ots::modifyBasisByOrientation(drv_div_basis,
1450  drv_div_basis_tmp,
1451  drv_orts,
1452  intrepid_basis);
1453  }
1454  if(build_weighted) {
1455  auto drv_div_basis = Kokkos::subview(weighted_div_basis.get_view(), range_cell, Kokkos::ALL(), Kokkos::ALL());
1456  auto drv_div_basis_tmp = Kokkos::createDynRankView(weighted_div_basis.get_view(),
1457  "drv_div_basis_tmp",
1458  drv_div_basis.extent(0), // C
1459  drv_div_basis.extent(1), // F
1460  drv_div_basis.extent(2)); // P
1461  Kokkos::deep_copy(drv_div_basis_tmp, drv_div_basis);
1462  ots::modifyBasisByOrientation(drv_div_basis,
1463  drv_div_basis_tmp,
1464  drv_orts,
1465  intrepid_basis);
1466  }
1467  }
1468  }
1469 }
1470 
1471 // method for applying orientations
1472 template <typename Scalar>
1474 applyOrientations(const PHX::MDField<const Scalar,Cell,BASIS> & orientations)
1475 {
1476 
1477  TEUCHOS_TEST_FOR_EXCEPT_MSG(true,"panzer::BasisValues2::applyOrientations : this should not be called.");
1478 
1479  int num_cell = orientations.extent(0);
1480  int num_basis = orientations.extent(1);
1481  int num_dim = basis_layout->dimension();
1482  int num_ip = basis_layout->numPoints();
1483  PureBasis::EElementSpace elmtspace = getElementSpace();
1484 
1485  if(elmtspace==PureBasis::HCURL && num_dim==2) {
1486 
1487  // setup the orientations for the trial space
1488  // Intrepid2::FunctionSpaceTools::applyFieldSigns<Scalar>(basis_vector,orientations);
1489 
1490  for (int c=0; c<num_cell; c++)
1491  for (int b=0; b<num_basis; b++)
1492  for (int p=0; p<num_ip; p++)
1493  for (int d=0; d<num_dim; d++)
1494  basis_vector(c, b, p, d) *= orientations(c, b);
1495 
1496  if(compute_derivatives) {
1497  // Intrepid2::FunctionSpaceTools::applyFieldSigns<Scalar>(curl_basis_scalar,orientations);
1498  for (int c=0; c<num_cell; c++)
1499  for (int b=0; b<num_basis; b++)
1500  for (int p=0; p<num_ip; p++)
1501  curl_basis_scalar(c, b, p) *= orientations(c, b);
1502  }
1503 
1504  // setup the orientations for the test space
1505  if(build_weighted) {
1506  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::applyFieldSigns(weighted_basis_vector.get_view(),orientations.get_view());
1507 
1508  if(compute_derivatives)
1509  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::applyFieldSigns(weighted_curl_basis_scalar.get_view(),orientations.get_view());
1510  }
1511  }
1512  else if(elmtspace==PureBasis::HCURL && num_dim==3) {
1513 
1514  // setup the orientations for the trial space
1515  // Intrepid2::FunctionSpaceTools::applyFieldSigns<Scalar>(basis_vector,orientations);
1516 
1517  for (int c=0; c<num_cell; c++)
1518  for (int b=0; b<num_basis; b++)
1519  for (int p=0; p<num_ip; p++)
1520  for (int d=0; d<num_dim; d++)
1521  basis_vector(c, b, p, d) *= orientations(c, b);
1522 
1523  if(compute_derivatives) {
1524  // Intrepid2::FunctionSpaceTools::applyFieldSigns<Scalar>(curl_basis_vector,orientations);
1525  for (int c=0; c<num_cell; c++)
1526  for (int b=0; b<num_basis; b++)
1527  for (int p=0; p<num_ip; p++)
1528  for (int d=0; d<num_dim; d++)
1529  curl_basis_vector(c, b, p,d) *= orientations(c, b);
1530  }
1531 
1532  // setup the orientations for the test space
1533  if(build_weighted) {
1534  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::applyFieldSigns(weighted_basis_vector.get_view(),orientations.get_view());
1535 
1536  if(compute_derivatives)
1537  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::applyFieldSigns(weighted_curl_basis_vector.get_view(),orientations.get_view());
1538  }
1539  }
1540  else if(elmtspace==PureBasis::HDIV) {
1541  // setup the orientations for the trial space
1542  // Intrepid2::FunctionSpaceTools::applyFieldSigns<Scalar>(basis_vector,orientations);
1543 
1544  for (int c=0; c<num_cell; c++)
1545  for (int b=0; b<num_basis; b++)
1546  for (int p=0; p<num_ip; p++)
1547  for (int d=0; d<num_dim; d++)
1548  basis_vector(c, b, p, d) *= orientations(c, b);
1549 
1550  if(compute_derivatives) {
1551  // Intrepid2::FunctionSpaceTools::applyFieldSigns<Scalar>(div_basis,orientations);
1552 
1553  for (int c=0; c<num_cell; c++)
1554  for (int b=0; b<num_basis; b++)
1555  for (int p=0; p<num_ip; p++)
1556  div_basis(c, b, p) *= orientations(c, b);
1557  }
1558 
1559  // setup the orientations for the test space
1560  if(build_weighted) {
1561  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::applyFieldSigns(weighted_basis_vector.get_view(),orientations.get_view());
1562 
1563  if(compute_derivatives)
1564  Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>::applyFieldSigns(weighted_div_basis.get_view(),orientations.get_view());
1565  }
1566  }
1567 }
1568 
1569 template <typename Scalar>
1571 { return basis_layout->getBasis()->getElementSpace(); }
1572 
1573 template <typename Scalar>
1576  bool computeDerivatives)
1577 {
1578  MDFieldArrayFactory af(prefix,alloc_arrays);
1579 
1580  compute_derivatives = computeDerivatives;
1581  basis_layout = layout;
1582  Teuchos::RCP<const panzer::PureBasis> basisDesc = layout->getBasis();
1583 
1584  // for convience pull out basis and quadrature information
1585  int num_quad = layout->numPoints();
1586  int dim = basisDesc->dimension();
1587  int card = basisDesc->cardinality();
1588  int numcells = basisDesc->numCells();
1589  panzer::PureBasis::EElementSpace elmtspace = basisDesc->getElementSpace();
1591 
1592  intrepid_basis = basisDesc->getIntrepid2Basis<PHX::Device::execution_space,Scalar,Scalar>();
1593 
1594  // allocate field containers
1595  // field sizes defined by http://trilinos.sandia.gov/packages/docs/dev/packages/intrepid/doc/html/basis_page.html#basis_md_array_sec
1596 
1597  // compute basis fields
1598  if(elmtspace==panzer::PureBasis::HGRAD) {
1599  // HGRAD is a nodal field
1600 
1601  // build values
1603  basis_ref_scalar = af.buildStaticArray<Scalar,BASIS,IP>("basis_ref",card,num_quad); // F, P
1604  basis_scalar = af.buildStaticArray<Scalar,Cell,BASIS,IP>("basis",numcells,card,num_quad);
1605 
1606  if(build_weighted)
1607  weighted_basis_scalar = af.buildStaticArray<Scalar,Cell,BASIS,IP>("weighted_basis",numcells,card,num_quad);
1608 
1609  // build gradients
1611 
1612  if(compute_derivatives) {
1613  grad_basis_ref = af.buildStaticArray<Scalar,BASIS,IP,Dim>("grad_basis_ref",card,num_quad,dim); // F, P, D
1614  grad_basis = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("grad_basis",numcells,card,num_quad,dim);
1615 
1616  if(build_weighted)
1617  weighted_grad_basis = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("weighted_grad_basis",numcells,card,num_quad,dim);
1618  }
1619 
1620  // build curl
1622 
1623  // nothing - HGRAD does not support CURL operation
1624  }
1625  else if(elmtspace==panzer::PureBasis::HCURL) {
1626  // HCURL is a vector field
1627 
1628  // build values
1630 
1631  basis_ref_vector = af.buildStaticArray<Scalar,BASIS,IP,Dim>("basis_ref",card,num_quad,dim); // F, P, D
1632  basis_vector = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("basis",numcells,card,num_quad,dim);
1633 
1634  if(build_weighted)
1635  weighted_basis_vector = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("weighted_basis",numcells,card,num_quad,dim);
1636 
1637  // build gradients
1639 
1640  // nothing - HCURL does not support GRAD operation
1641 
1642  // build curl
1644 
1645  if(compute_derivatives) {
1646  if(dim==2) {
1647  // curl of HCURL basis is not dimension dependent
1648  curl_basis_ref_scalar = af.buildStaticArray<Scalar,BASIS,IP>("curl_basis_ref",card,num_quad); // F, P
1649  curl_basis_scalar = af.buildStaticArray<Scalar,Cell,BASIS,IP>("curl_basis",numcells,card,num_quad);
1650 
1651  if(build_weighted)
1652  weighted_curl_basis_scalar = af.buildStaticArray<Scalar,Cell,BASIS,IP>("weighted_curl_basis",numcells,card,num_quad);
1653  }
1654  else if(dim==3){
1655  curl_basis_ref_vector = af.buildStaticArray<Scalar,BASIS,IP,Dim>("curl_basis_ref",card,num_quad,dim); // F, P, D
1656  curl_basis_vector = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("curl_basis",numcells,card,num_quad,dim);
1657 
1658  if(build_weighted)
1659  weighted_curl_basis_vector = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("weighted_curl_basis",numcells,card,num_quad,dim);
1660  }
1661  else { TEUCHOS_ASSERT(false); } // what do I do with 1D?
1662  }
1663  }
1664  else if(elmtspace==panzer::PureBasis::HDIV) {
1665  // HDIV is a vector field
1666 
1667  // build values
1669 
1670  basis_ref_vector = af.buildStaticArray<Scalar,BASIS,IP,Dim>("basis_ref",card,num_quad,dim); // F, P, D
1671  basis_vector = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("basis",numcells,card,num_quad,dim);
1672 
1673  if(build_weighted)
1674  weighted_basis_vector = af.buildStaticArray<Scalar,Cell,BASIS,IP,Dim>("weighted_basis",numcells,card,num_quad,dim);
1675 
1676  // build gradients
1678 
1679  // nothing - HCURL does not support GRAD operation
1680 
1681  // build curl
1683 
1684  // nothing - HDIV does not support CURL operation
1685 
1686  // build div
1688 
1689  if(compute_derivatives) {
1690  div_basis_ref = af.buildStaticArray<Scalar,BASIS,IP>("div_basis_ref",card,num_quad); // F, P
1691  div_basis = af.buildStaticArray<Scalar,Cell,BASIS,IP>("div_basis",numcells,card,num_quad);
1692 
1693  if(build_weighted)
1694  weighted_div_basis = af.buildStaticArray<Scalar,Cell,BASIS,IP>("weighted_div_basis",numcells,card,num_quad);
1695  }
1696  }
1697  else if(elmtspace==panzer::PureBasis::CONST || elmtspace==panzer::PureBasis::HVOL) {
1698  // CONST is a nodal field
1699 
1700  // build values
1702  basis_ref_scalar = af.buildStaticArray<Scalar,BASIS,IP>("basis_ref",card,num_quad); // F, P
1703  basis_scalar = af.buildStaticArray<Scalar,Cell,BASIS,IP>("basis",numcells,card,num_quad);
1704 
1705  if(build_weighted)
1706  weighted_basis_scalar = af.buildStaticArray<Scalar,Cell,BASIS,IP>("weighted_basis",numcells,card,num_quad);
1707 
1708  // build gradients
1710 
1711  // nothing - CONST does not support GRAD operation
1712 
1713  // build curl
1715 
1716  // nothing - CONST does not support CURL operation
1717 
1718  // build div
1720 
1721  // nothing - CONST does not support DIV operation
1722  }
1723  else { TEUCHOS_ASSERT(false); }
1724 
1725  basis_coordinates_ref = af.buildStaticArray<Scalar,BASIS,Dim>("basis_coordinates_ref",card,dim);
1726  basis_coordinates = af.buildStaticArray<Scalar,Cell,BASIS,Dim>("basis_coordinates",numcells,card,dim);
1727 }
1728 
1729 } // namespace panzer
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack...dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
void evaluateReferenceValues(const PHX::MDField< Scalar, IP, Dim > &cub_points, bool compute_derivatives, bool use_vertex_coordinates)
EElementSpace getElementSpace() const
int cardinality() const
Returns the number of basis coefficients.
ArrayTraits< Scalar, PHX::MDField< Scalar, void, void, void, void, void, void, void, void > >::size_type size_type
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
Teuchos::RCP< const PureBasis > getBasis() const
int dimension() const
Returns the dimension of the basis from the topology.
void applyOrientations(const PHX::MDField< const Scalar, Cell, BASIS > &orientations)
Method to apply orientations to a basis values container.
Teuchos::RCP< const shards::CellTopology > getCellTopology() const
#define TEUCHOS_TEST_FOR_EXCEPT_MSG(throw_exception_test, msg)
void evaluateValues(const PHX::MDField< Scalar, IP, Dim, void, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac, const PHX::MDField< Scalar, Cell, IP, void, void, void, void, void, void > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv, const int in_num_cells=-1)
void evaluateBasisCoordinates(const PHX::MDField< Scalar, Cell, NODE, Dim > &vertex_coordinates, const int in_num_cells=-1)
void evaluateValues_HVol(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, void, void, void, void, void, void > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv, const PHX::MDField< Scalar, Cell, IP > &weighted_measure, const int in_num_cells)
void evaluateValues_HGrad(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv, const PHX::MDField< Scalar, Cell, IP > &weighted_measure, const int in_num_cells)
void evaluateValues_Const(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv, const PHX::MDField< Scalar, Cell, IP > &weighted_measure, const int in_num_cells)
PureBasis::EElementSpace getElementSpace() const
Teuchos::RCP< Intrepid2::Basis< PHX::Device::execution_space, double, double > > getIntrepid2Basis() const
void evaluateValues_HDiv(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac, const PHX::MDField< Scalar, Cell, IP, void, void, void, void, void, void > &jac_det, const PHX::MDField< Scalar, Cell, IP > &weighted_measure, const int in_num_cells)
void setupArrays(const Teuchos::RCP< const panzer::BasisIRLayout > &basis, bool computeDerivatives=true)
Sizes/allocates memory for arrays.
int numCells() const
Returns the number of cells in the data layouts.
Description and data layouts associated with a particular basis.
#define TEUCHOS_ASSERT(assertion_test)
void evaluateValuesCV(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cell_cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac, const PHX::MDField< Scalar, Cell, IP, void, void, void, void, void, void > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv)
PHX::MDField< Scalar > buildArray(const std::string &str, int d0) const
void evaluateValues_HCurl(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac, const PHX::MDField< Scalar, Cell, IP, void, void, void, void, void, void > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv, const PHX::MDField< Scalar, Cell, IP > &weighted_measure, const int in_num_cells)
PHX::MDField< Scalar > ArrayDynamic