Intrepid2
Intrepid2_FunctionSpaceToolsDef.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
16 #ifndef __INTREPID2_FUNCTIONSPACETOOLS_DEF_HPP__
17 #define __INTREPID2_FUNCTIONSPACETOOLS_DEF_HPP__
18 
21 
22 #include "Teuchos_TimeMonitor.hpp"
23 
24 namespace Intrepid2 {
25  // ------------------------------------------------------------------------------------
26  template<typename DeviceType>
27  template<typename outputValueType, class ...outputProperties,
28  typename inputValueType, class ...inputProperties>
29  void
31  HGRADtransformVALUE( Kokkos::DynRankView<outputValueType,outputProperties...> output,
32  const Kokkos::DynRankView<inputValueType, inputProperties...> input ) {
33  if(output.rank() == input.rank()) {
34 #ifdef HAVE_INTREPID2_DEBUG
35  {
36  for (size_type i=0;i< input.rank();++i) {
37  INTREPID2_TEST_FOR_EXCEPTION( (input.extent(i) != output.extent(i)), std::invalid_argument,
38  ">>> ERROR (FunctionSpaceTools::HGRADtransformVALUE): Dimensions of input and output fields containers do not match.");
39  }
40  }
41 #endif
42  RealSpaceTools<DeviceType>::clone(output, input);
43  }
44  else
46  }
47 
48  template<typename DeviceType>
49  template<typename outputValueType, class ...outputProperties,
50  typename inputValueType, class ...inputProperties>
51  void
53  mapHGradDataFromPhysToRef( Kokkos::DynRankView<outputValueType,outputProperties...> output,
54  const Kokkos::DynRankView<inputValueType, inputProperties...> input ) {
55  if(output.rank() == input.rank()) {
56 #ifdef HAVE_INTREPID2_DEBUG
57  {
58  for (size_type i=0;i< input.rank();++i) {
59  INTREPID2_TEST_FOR_EXCEPTION( (input.extent(i) != output.extent(i)), std::invalid_argument,
60  ">>> ERROR (FunctionSpaceTools::mapHGradDataFromPhysToRef): Dimensions of input and output fields containers do not match.");
61  }
62  }
63 #endif
64  RealSpaceTools<DeviceType>::clone(output, input);
65  }
66  else
68  }
69 
70  template<typename DeviceType>
71  template<typename outputValueType, class ...outputProperties,
72  typename inputValueType, class ...inputProperties>
73  void
75  mapHGradDataFromPhysSideToRefSide( Kokkos::DynRankView<outputValueType,outputProperties...> output,
76  const Kokkos::DynRankView<inputValueType, inputProperties...> input ) {
77  if(output.rank() == input.rank()) {
78 #ifdef HAVE_INTREPID2_DEBUG
79  {
80  for (size_type i=0;i< input.rank();++i) {
81  INTREPID2_TEST_FOR_EXCEPTION( (input.extent(i) != output.extent(i)), std::invalid_argument,
82  ">>> ERROR (FunctionSpaceTools::mapHGradDataSideFromPhysToRefSide): Dimensions of input and output fields containers do not match.");
83  }
84  }
85 #endif
86  RealSpaceTools<DeviceType>::clone(output, input);
87  }
88  else
90  }
91 
92  // ------------------------------------------------------------------------------------
93 
94  namespace FunctorFunctionSpaceTools {
98  }
99 
100  template<typename DeviceType>
101  template<typename OutputValViewType,
102  typename JacobianInverseViewType,
103  typename InputValViewType>
104  void
106  HGRADtransformGRAD( OutputValViewType outputVals,
107  const JacobianInverseViewType jacobianInverse,
108  const InputValViewType inputVals ) {
109  return HCURLtransformVALUE(outputVals, jacobianInverse, inputVals);
110  }
111 
112  // ------------------------------------------------------------------------------------
113 
114  template<typename DeviceType>
115  template<typename outputValValueType, class ...outputValProperties,
116  typename jacobianInverseValueType, class ...jacobianInverseProperties,
117  typename inputValValueType, class ...inputValProperties>
118  void
120  HCURLtransformVALUE( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
121  const Kokkos::DynRankView<jacobianInverseValueType,jacobianInverseProperties...> jacobianInverse,
122  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
123  ArrayTools<DeviceType>::matvecProductDataField(outputVals, jacobianInverse, inputVals, 'T');
124  }
125 
126  template<typename DeviceType>
127  template<typename outputValValueType, class ...outputValProperties,
128  typename jacobianValueType, class ...jacobianProperties,
129  typename inputValValueType, class ...inputValProperties>
130  void
132  mapHCurlDataFromPhysToRef( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
133  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
134  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
135  ArrayTools<DeviceType>::matvecProductDataData(outputVals, jacobian, inputVals, 'T');
136  }
137 
138 
139  namespace FunctorFunctionSpaceTools {
140 
141  template<typename outViewType,
142  typename inputViewType,
143  typename metricViewType
144  >
146  outViewType output_;
147  const inputViewType input_;
148  const metricViewType metricTensorDet_;
149 
150  KOKKOS_INLINE_FUNCTION
151  F_negativeWeighted2dInputCrossK( outViewType output,
152  const inputViewType input,
153  const metricViewType metricTensorDet)
154  : output_(output), input_(input), metricTensorDet_(metricTensorDet){};
155 
156  KOKKOS_INLINE_FUNCTION
157  void operator()(const size_type ic) const {
158  for (size_t pt=0; pt < input_.extent(1); pt++) {
159  auto measure = std::sqrt(metricTensorDet_(ic,pt));
160  output_(ic,pt,0) = - measure*input_(ic,pt,1);
161  output_(ic,pt,1) = measure*input_(ic,pt,0);
162  }
163  }
164  };
165  }
166 
167  template<typename DeviceType>
168  template<typename outputValValueType, class ...outputValProperties,
169  typename tangentsValueType, class ...tangentsProperties,
170  typename metricTensorInvValueType, class ...metricTensorInvProperties,
171  typename metricTensorDetValueType, class ...metricTensorDetProperties,
172  typename inputValValueType, class ...inputValProperties>
173  void
175  mapHCurlDataCrossNormalFromPhysSideToRefSide( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
176  const Kokkos::DynRankView<tangentsValueType, tangentsProperties...> tangents,
177  const Kokkos::DynRankView<metricTensorInvValueType,metricTensorInvProperties...> metricTensorInv,
178  const Kokkos::DynRankView<metricTensorDetValueType,metricTensorDetProperties...> metricTensorDet,
179  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
180  auto work = Kokkos::DynRankView<outputValValueType, outputValProperties...>("work", outputVals.extent(0), outputVals.extent(1),outputVals.extent(2));
181  ArrayTools<DeviceType>::matvecProductDataData(outputVals, tangents, inputVals, 'T');
182  typename DeviceType::execution_space().fence();
183  ArrayTools<DeviceType>::matvecProductDataData(work, metricTensorInv, outputVals);
184  typename DeviceType::execution_space().fence();
186  Kokkos::parallel_for(Kokkos::RangePolicy<ExecSpaceType>(0,outputVals.extent(0)), FunctorType(outputVals, work, metricTensorDet) );
187  }
188 
189 
190  namespace FunctorFunctionSpaceTools {
191 
192  template<typename outViewType,
193  typename inputViewType,
194  typename metricViewType
195  >
196  struct F_weighedInput {
197  outViewType output_;
198  const inputViewType input_;
199  const metricViewType metricTensorDet_;
200  const double scaling_;
201 
202  KOKKOS_INLINE_FUNCTION
203  F_weighedInput( outViewType output,
204  const inputViewType input,
205  const metricViewType metricTensorDet,
206  const double scaling)
207  : output_(output), input_(input), metricTensorDet_(metricTensorDet), scaling_(scaling){};
208 
209  KOKKOS_INLINE_FUNCTION
210  void operator()(const size_type ic) const {
211  for (size_t pt=0; pt < input_.extent(1); pt++) {
212  auto measure = std::sqrt(metricTensorDet_(ic,pt));
213  output_.access(ic,pt) = scaling_ * measure * input_.access(ic,pt);
214  }
215  }
216  };
217  }
218 
219  template<typename DeviceType>
220  template<typename outputValValueType, class ...outputValProperties,
221  typename jacobianDetValueType, class ...jacobianDetProperties,
222  typename inputValValueType, class ...inputValProperties>
223  void
226  Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
227  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> metricTensorDet,
228  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
230  Kokkos::parallel_for(Kokkos::RangePolicy<ExecSpaceType>(0,outputVals.extent(0)), FunctorType(outputVals, inputVals, metricTensorDet, -1.0) );
231  }
232 
233  // ------------------------------------------------------------------------------------
234 
235  template<typename DeviceType>
236  template<typename outputValValueType, class ...outputValProperties,
237  typename jacobianValueType, class ...jacobianProperties,
238  typename jacobianDetValueType, class ...jacobianDetProperties,
239  typename inputValValueType, class ...inputValProperties>
240  void
242  HCURLtransformCURL( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
243  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
244  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
245  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
246  if(jacobian.data()==NULL || jacobian.extent(2)==2) //2D case
247  return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
248  else
249  return HDIVtransformVALUE(outputVals, jacobian, jacobianDet, inputVals);
250  }
251 
252  // ------------------------------------------------------------------------------------
253 
254  template<typename DeviceType>
255  template<typename outputValValueType, class ...outputValProperties,
256  typename jacobianDetValueType, class ...jacobianDetProperties,
257  typename inputValValueType, class ...inputValProperties>
258  void
260  HCURLtransformCURL( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
261  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
262  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
263 #ifdef HAVE_INTREPID2_DEBUG
264  {
265  INTREPID2_TEST_FOR_EXCEPTION( outputVals.rank() == 4, std::invalid_argument,
266  ">>> ERROR (FunctionSpaceTools::HCURLtransformCURL): Output rank must have rank 3.\n If these are 3D fields, then use the appropriate overload of this function.");
267  }
268 #endif
269  return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
270  }
271 
272  // ------------------------------------------------------------------------------------
273 
274  template<typename DeviceType>
275  template<typename outputValValueType, class ...outputValProperties,
276  typename jacobianValueType, class ...jacobianProperties,
277  typename jacobianDetValueType, class ...jacobianDetProperties,
278  typename inputValValueType, class ...inputValProperties>
279  void
281  HGRADtransformCURL( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
282  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
283  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
284  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
285 #ifdef HAVE_INTREPID2_DEBUG
286  {
287  INTREPID2_TEST_FOR_EXCEPTION( outputVals.extent(3)!=2, std::invalid_argument,
288  ">>> ERROR (FunctionSpaceTools::HGRADtransformCURL):\n output field is 3D by the function is meant for 2D fields");
289  }
290 #endif
291  return HDIVtransformVALUE(outputVals, jacobian, jacobianDet, inputVals);
292  }
293 
294  // ------------------------------------------------------------------------------------
295 
296  template<typename DeviceType>
297  template<typename outputValValueType, class ...outputValProperties,
298  typename jacobianValueType, class ...jacobianProperties,
299  typename jacobianDetValueType, class ...jacobianDetProperties,
300  typename inputValValueType, class ...inputValProperties>
301  void
303  HDIVtransformVALUE( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
304  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
305  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
306  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
307  ArrayTools<DeviceType>::matvecProductDataField(outputVals, jacobian, inputVals, 'N');
308  ArrayTools<DeviceType>::scalarMultiplyDataField(outputVals, jacobianDet, outputVals, true);
309  }
310 
311  template<typename DeviceType>
312  template<typename outputValValueType, class ...outputValProperties,
313  typename jacobianInverseValueType, class ...jacobianInverseProperties,
314  typename jacobianDetValueType, class ...jacobianDetProperties,
315  typename inputValValueType, class ...inputValProperties>
316  void
318  mapHDivDataFromPhysToRef( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
319  const Kokkos::DynRankView<jacobianInverseValueType, jacobianInverseProperties...> jacobianInv,
320  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
321  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
322  ArrayTools<DeviceType>::matvecProductDataData(outputVals, jacobianInv, inputVals, 'N');
323  typename DeviceType::execution_space().fence();
324  ArrayTools<DeviceType>::scalarMultiplyDataData(outputVals, jacobianDet, outputVals, false);
325  }
326 
327 
328 
329  template<typename DeviceType>
330  template<typename outputValValueType, class ...outputValProperties,
331  typename jacobianDetValueType, class ...jacobianDetProperties,
332  typename inputValValueType, class ...inputValProperties>
333  void
336  Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
337  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> metricTensorDet,
338  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
340  Kokkos::parallel_for(Kokkos::RangePolicy<ExecSpaceType>(0,outputVals.extent(0)), FunctorType(outputVals, inputVals, metricTensorDet, 1.0) );
341  }
342 
343  // ------------------------------------------------------------------------------------
344 
345  template<typename DeviceType>
346  template<typename outputValValueType, class ...outputValProperties,
347  typename jacobianDetValueType, class ...jacobianDetProperties,
348  typename inputValValueType, class ...inputValProperties>
349  void
351  HDIVtransformDIV( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
352  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
353  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
354  return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
355  }
356 
357  // ------------------------------------------------------------------------------------
358 
359  template<typename DeviceType>
360  template<typename outputValValueType, class ...outputValProperties,
361  typename jacobianDetValueType, class ...jacobianDetProperties,
362  typename inputValValueType, class ...inputValProperties>
363  void
365  HVOLtransformVALUE( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
366  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
367  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
368  ArrayTools<DeviceType>::scalarMultiplyDataField(outputVals, jacobianDet, inputVals, true);
369  }
370 
371  template<typename DeviceType>
372  template<typename outputValValueType, class ...outputValProperties,
373  typename jacobianDetValueType, class ...jacobianDetProperties,
374  typename inputValValueType, class ...inputValProperties>
375  void
377  mapHVolDataFromPhysToRef( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
378  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
379  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
380  ArrayTools<DeviceType>::scalarMultiplyDataData(outputVals, jacobianDet, inputVals, false);
381  }
382 
383 
384 
385  // ------------------------------------------------------------------------------------
386 
387  template<typename DeviceType>
388  template<typename outputValueValueType, class ...outputValueProperties,
389  typename leftValueValueType, class ...leftValueProperties,
390  typename rightValueValueType, class ...rightValueProperties>
391  void
393  integrate( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
394  const Kokkos::DynRankView<leftValueValueType, leftValueProperties...> leftValues,
395  const Kokkos::DynRankView<rightValueValueType, rightValueProperties...> rightValues,
396  const bool sumInto ) {
397 
398 #ifdef HAVE_INTREPID2_DEBUG
399  {
400  INTREPID2_TEST_FOR_EXCEPTION( leftValues.rank() < 2 ||
401  leftValues.rank() > 4, std::invalid_argument,
402  ">>> ERROR (FunctionSpaceTools::integrate): Left data must have rank 2, 3 or 4.");
403  INTREPID2_TEST_FOR_EXCEPTION( outputValues.rank() < 1 ||
404  outputValues.rank() > 3, std::invalid_argument,
405  ">>> ERROR (FunctionSpaceTools::integrate): Output values must have rank 1, 2 or 3.");
406  }
407 #endif
408 
409  const ordinal_type outRank = outputValues.rank();
410  const ordinal_type leftRank = leftValues.rank();
411  const ordinal_type mode = outRank*10 + leftRank;
412 
413  switch (mode) {
414  // DataData
415  case 12:
417  leftValues,
418  rightValues,
419  sumInto );
420  break;
421  case 13:
423  leftValues,
424  rightValues,
425  sumInto );
426  break;
427  case 14:
429  leftValues,
430  rightValues,
431  sumInto );
432  break;
433 
434  // DataField
435  case 22:
437  leftValues,
438  rightValues,
439  sumInto );
440  break;
441  case 23:
443  leftValues,
444  rightValues,
445  sumInto );
446  break;
447  case 24:
449  leftValues,
450  rightValues,
451  sumInto );
452  break;
453 
454  // FieldField
455  case 33:
457  leftValues,
458  rightValues,
459  sumInto );
460  break;
461  case 34:
463  leftValues,
464  rightValues,
465  sumInto );
466  break;
467  case 35:
469  leftValues,
470  rightValues,
471  sumInto );
472  break;
473  default: {
474  INTREPID2_TEST_FOR_EXCEPTION( outRank < 1 || outRank > 3, std::runtime_error,
475  ">>> ERROR (FunctionSpaceTools::integrate): outRank must be 1,2, or 3.");
476  INTREPID2_TEST_FOR_EXCEPTION( leftRank < 2 || leftRank > 4, std::runtime_error,
477  ">>> ERROR (FunctionSpaceTools::integrate): leftRank must be 1,2, 3 or 4.");
478  }
479  }
480  }
481 
482  // ------------------------------------------------------------------------------------
483  namespace FunctorFunctionSpaceTools {
487  template<typename outputValViewType,
488  typename inputDetViewType,
489  typename inputWeightViewType>
491  outputValViewType _outputVals;
492  const inputDetViewType _inputDet;
493  const inputWeightViewType _inputWeight;
494 
495  KOKKOS_INLINE_FUNCTION
496  F_computeCellMeasure( outputValViewType outputVals_,
497  inputDetViewType inputDet_,
498  inputWeightViewType inputWeight_)
499  : _outputVals(outputVals_),
500  _inputDet(inputDet_),
501  _inputWeight(inputWeight_) {}
502 
503  typedef ordinal_type value_type;
504 
505 // KOKKOS_INLINE_FUNCTION
506 // void init( value_type &dst ) const {
507 // dst = false;
508 // }
509 
510 // KOKKOS_INLINE_FUNCTION
511 // void join( volatile value_type &dst,
512 // const volatile value_type &src ) const {
513 // dst |= src;
514 // }
515 
516  KOKKOS_INLINE_FUNCTION
517  void operator()(const size_type cl, value_type &dst) const {
518  // negative jacobian check
519  const bool hasNegativeDet = (_inputDet(cl, 0) < 0.0);
520  dst |= hasNegativeDet;
521 
522  // make it absolute
523  const auto sign = (hasNegativeDet ? -1.0 : 1.0);
524  const ordinal_type pt_end = _outputVals.extent(1);
525  for (ordinal_type pt=0;pt<pt_end;++pt)
526  _outputVals(cl, pt) = sign*_inputDet(cl, pt)*_inputWeight(pt);
527  }
528  };
529  }
530 
531  template<typename DeviceType>
532  template<typename OutputValViewType,
533  typename InputDetViewType,
534  typename InputWeightViewType>
535  bool
537  computeCellMeasure( OutputValViewType outputVals,
538  const InputDetViewType inputDet,
539  const InputWeightViewType inputWeights ) {
540 #ifdef HAVE_INTREPID2_DEBUG
541  {
542  INTREPID2_TEST_FOR_EXCEPTION( rank(inputDet) != 2 ||
543  rank(inputWeights) != 1 ||
544  rank(outputVals) != 2, std::invalid_argument,
545  ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Ranks are not compatible.");
546  INTREPID2_TEST_FOR_EXCEPTION( outputVals.extent(0) != inputDet.extent(0), std::invalid_argument,
547  ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Cell dimension does not match.");
548  INTREPID2_TEST_FOR_EXCEPTION( inputDet.extent(1) != outputVals.extent(1) ||
549  inputWeights.extent(0) != outputVals.extent(1), std::invalid_argument,
550  ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Point dimension does not match.");
551  }
552 #endif
553  constexpr bool are_accessible =
554  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
555  typename decltype(outputVals)::memory_space>::accessible &&
556  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
557  typename decltype(inputDet)::memory_space>::accessible &&
558  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
559  typename decltype(inputWeights)::memory_space>::accessible;
560  static_assert(are_accessible, "FunctionSpaceTools<DeviceType>::computeCellMeasure(..): input/output views' memory spaces are not compatible with DeviceType");
561 
563  <decltype(outputVals),decltype(inputDet),decltype(inputWeights)>;
564 
565  const ordinal_type C = inputDet.extent(0);
566  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
567 
568  typename FunctorType::value_type hasNegativeDet = false;
569  Kokkos::parallel_reduce( policy, FunctorType(outputVals, inputDet, inputWeights), hasNegativeDet );
570 
571  return hasNegativeDet;
572  }
573 
574  // ------------------------------------------------------------------------------------
575 
576  template<typename DeviceType>
577  template<typename outputValValueType, class ...outputValProperties,
578  typename inputJacValueType, class ...inputJacProperties,
579  typename inputWeightValueType, class ...inputWeightProperties,
580  typename scratchValueType, class ...scratchProperties>
581  void
583  computeFaceMeasure( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
584  const Kokkos::DynRankView<inputJacValueType, inputJacProperties...> inputJac,
585  const Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeights,
586  const ordinal_type whichFace,
587  const shards::CellTopology parentCell,
588  const Kokkos::DynRankView<scratchValueType, scratchProperties...> scratch ) {
589 #ifdef HAVE_INTREPID2_DEBUG
590  INTREPID2_TEST_FOR_EXCEPTION( inputJac.rank() != 4, std::invalid_argument,
591  ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Input Jacobian container must have rank 4.");
592  INTREPID2_TEST_FOR_EXCEPTION( scratch.rank() != 1, std::invalid_argument,
593  ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Scratch view imust have rank 1.");
594  INTREPID2_TEST_FOR_EXCEPTION( scratch.size() < inputJac.size(), std::invalid_argument,
595  ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Scratch storage must be greater than or equal to inputJac's one.");
596 #endif
597 
598  // face normals (reshape scratch)
599  // Kokkos::DynRankView<scratchValueType,scratchProperties...> faceNormals(scratch.data(),
600  // inputJac.extent(0),
601  // inputJac.extent(1),
602  // inputJac.extent(2));
603  auto vcprop = Kokkos::common_view_alloc_prop(scratch);
604  //typedef Kokkos::DynRankView<scratchValueType, typename decltype(scratch)::memory_space> viewType;
605  typedef Kokkos::DynRankView<scratchValueType, DeviceType> viewType;
606  viewType faceNormals(Kokkos::view_wrap(scratch.data(), vcprop),
607  inputJac.extent(0),
608  inputJac.extent(1),
609  inputJac.extent(2));
610 
611  // compute normals
612  CellTools<DeviceType>::getPhysicalFaceNormals(faceNormals, inputJac, whichFace, parentCell);
613 
614  // compute lenghts of normals
615  RealSpaceTools<DeviceType>::vectorNorm(outputVals, faceNormals, NORM_TWO);
616 
617  // multiply with weights
618  ArrayTools<DeviceType>::scalarMultiplyDataData(outputVals, outputVals, inputWeights);
619  }
620 
621  // ------------------------------------------------------------------------------------
622 
623  template<typename DeviceType>
624  template<typename outputValValueType, class ...outputValProperties,
625  typename inputJacValueType, class ...inputJacProperties,
626  typename inputWeightValueType, class ...inputWeightProperties,
627  typename scratchValueType, class ...scratchProperties>
628  void
630  computeEdgeMeasure( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
631  const Kokkos::DynRankView<inputJacValueType, inputJacProperties...> inputJac,
632  const Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeights,
633  const ordinal_type whichEdge,
634  const shards::CellTopology parentCell,
635  const Kokkos::DynRankView<scratchValueType, scratchProperties...> scratch ) {
636 #ifdef HAVE_INTREPID2_DEBUG
637  INTREPID2_TEST_FOR_EXCEPTION( (inputJac.rank() != 4), std::invalid_argument,
638  ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Input Jacobian container must have rank 4.");
639  INTREPID2_TEST_FOR_EXCEPTION( scratch.rank() != 1, std::invalid_argument,
640  ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Scratch view must have a rank 1.");
641  INTREPID2_TEST_FOR_EXCEPTION( scratch.size() < inputJac.size(), std::invalid_argument,
642  ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Scratch storage must be greater than or equal to inputJac'one.");
643 #endif
644 
645  // edge tangents (reshape scratch)
646  // Kokkos::DynRankView<scratchValueType,scratchProperties...> edgeTangents(scratch.data(),
647  // inputJac.extent(0),
648  // inputJac.extent(1),
649  // inputJac.extent(2));
650  auto vcprop = Kokkos::common_view_alloc_prop(scratch);
651  //typedef Kokkos::DynRankView<scratchValueType, typename decltype(scratch)::memory_space> viewType;
652  typedef Kokkos::DynRankView<scratchValueType, DeviceType> viewType;
653  viewType edgeTangents(Kokkos::view_wrap(scratch.data(), vcprop),
654  inputJac.extent(0),
655  inputJac.extent(1),
656  inputJac.extent(2));
657 
658  // compute normals
659  CellTools<DeviceType>::getPhysicalEdgeTangents(edgeTangents, inputJac, whichEdge, parentCell);
660 
661  // compute lenghts of tangents
662  RealSpaceTools<DeviceType>::vectorNorm(outputVals, edgeTangents, NORM_TWO);
663 
664  // multiply with weights
665  ArrayTools<DeviceType>::scalarMultiplyDataData(outputVals, outputVals, inputWeights);
666  }
667 
668  // ------------------------------------------------------------------------------------
669 
670  template<typename DeviceType>
671  template<typename outputValValueType, class ...outputValProperties,
672  typename inputMeasureValueType, class ...inputMeasureProperties,
673  typename inputValValueType, class ...inputValProperties>
674  void
676  multiplyMeasure( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
677  const Kokkos::DynRankView<inputMeasureValueType,inputMeasureProperties...> inputMeasure,
678  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
680  inputMeasure,
681  inputVals,
682  false );
683  }
684 
685  // ------------------------------------------------------------------------------------
686 
687  template<typename DeviceType>
688  template<typename outputFieldValueType, class ...outputFieldProperties,
689  typename inputDataValueType, class ...inputDataProperties,
690  typename inputFieldValueType, class ...inputFieldProperties>
691  void
693  scalarMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
694  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
695  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields,
696  const bool reciprocal ) {
698  inputData,
699  inputFields,
700  reciprocal );
701  }
702 
703  // ------------------------------------------------------------------------------------
704 
705  template<typename DeviceType>
706  template<typename outputDataValuetype, class ...outputDataProperties,
707  typename inputDataLeftValueType, class ...inputDataLeftProperties,
708  typename inputDataRightValueType, class ...inputDataRightProperties>
709  void
711  scalarMultiplyDataData( Kokkos::DynRankView<outputDataValuetype, outputDataProperties...> outputData,
712  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
713  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight,
714  const bool reciprocal ) {
716  inputDataLeft,
717  inputDataRight,
718  reciprocal );
719  }
720 
721  // ------------------------------------------------------------------------------------
722 
723  template<typename DeviceType>
724  template<typename outputFieldValueType, class ...outputFieldProperties,
725  typename inputDataValueType, class ...inputDataProperties,
726  typename inputFieldValueType, class ...inputFieldProperties>
727  void
729  dotMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
730  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
731  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields ) {
733  inputData,
734  inputFields );
735  }
736 
737  // ------------------------------------------------------------------------------------
738 
739  template<typename DeviceType>
740  template<typename outputDataValueType, class ...outputDataProperties,
741  typename inputDataLeftValueType, class ...inputDataLeftProperties,
742  typename inputDataRightValueType, class ...inputDataRightProperties>
743  void
745  dotMultiplyDataData( Kokkos::DynRankView<outputDataValueType, outputDataProperties...> outputData,
746  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
747  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight ) {
749  inputDataLeft,
750  inputDataRight );
751  }
752 
753  // ------------------------------------------------------------------------------------
754 
755  template<typename DeviceType>
756  template<typename outputFieldValueType, class ...outputFieldProperties,
757  typename inputDataValueType, class ...inputDataProperties,
758  typename inputFieldValueType, class ...inputFieldProperties>
759  void
761  vectorMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
762  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
763  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields ) {
764  const auto outRank = outputFields.rank();
765  switch (outRank) {
766  case 3:
767  case 4:
769  inputData,
770  inputFields );
771  break;
772  case 5:
774  inputData,
775  inputFields );
776  break;
777  default: {
778  INTREPID2_TEST_FOR_EXCEPTION( outRank < 3 && outRank > 5, std::runtime_error,
779  ">>> ERROR (FunctionSpaceTools::vectorMultiplyDataField): Output container must have rank 3, 4 or 5.");
780  }
781  }
782  }
783 
784  // ------------------------------------------------------------------------------------
785 
786  template<typename DeviceType>
787  template<typename outputDataValueType, class ...outputDataProperties,
788  typename inputDataLeftValueType, class ...inputDataLeftProperties,
789  typename inputDataRightValueType, class ...inputDataRightProperties>
790  void
792  vectorMultiplyDataData( Kokkos::DynRankView<outputDataValueType, outputDataProperties...> outputData,
793  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
794  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight ) {
795  const auto outRank = outputData.rank();
796  switch (outRank) {
797  case 2:
798  case 3:
800  inputDataLeft,
801  inputDataRight );
802  break;
803  case 4:
805  inputDataLeft,
806  inputDataRight );
807  break;
808  default: {
809  INTREPID2_TEST_FOR_EXCEPTION( outRank < 2 && outRank > 4, std::runtime_error,
810  ">>> ERROR (FunctionSpaceTools::vectorMultiplyDataData): Output container must have rank 2, 3 or 4.");
811  }
812  }
813  }
814 
815  // ------------------------------------------------------------------------------------
816 
817  template<typename DeviceType>
818  template<typename outputFieldValueType, class ...outputFieldProperties,
819  typename inputDataValueType, class ...inputDataProperties,
820  typename inputFieldValueType, class ...inputFieldProperties>
821  void
823  tensorMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
824  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
825  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields,
826  const char transpose ) {
827 
828  const auto outRank = outputFields.rank();
829  switch (outRank) {
830  case 4:
832  inputData,
833  inputFields,
834  transpose );
835  break;
836  case 5:
838  inputData,
839  inputFields,
840  transpose );
841  break;
842  default: {
843  INTREPID2_TEST_FOR_EXCEPTION( outRank < 4 && outRank > 5, std::runtime_error,
844  ">>> ERROR (FunctionSpaceTools::tensorMultiplyDataField): Output container must have rank 4 or 5.");
845  }
846  }
847  }
848 
849  // ------------------------------------------------------------------------------------
850 
851  template<typename DeviceType>
852  template<typename outputDataValueType, class ...outputDataProperties,
853  typename inputDataLeftValueType, class ...inputDataLeftProperties,
854  typename inputDataRightValueType, class ...inputDataRightProperties>
855  void
857  tensorMultiplyDataData( Kokkos::DynRankView<outputDataValueType, outputDataProperties...> outputData,
858  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
859  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight,
860  const char transpose ) {
861  const auto outRank = outputData.rank();
862  switch (outRank) {
863  case 3:
865  inputDataLeft,
866  inputDataRight,
867  transpose );
868  break;
869  case 4:
871  inputDataLeft,
872  inputDataRight,
873  transpose );
874  break;
875  default: {
876  INTREPID2_TEST_FOR_EXCEPTION( outRank < 4 && outRank > 5, std::runtime_error,
877  ">>> ERROR (FunctionSpaceTools::tensorMultiplyDataField): Output container must have rank 4 or 5.");
878  }
879  }
880  }
881 
882  // ------------------------------------------------------------------------------------
883 
884  namespace FunctorFunctionSpaceTools {
885 
889  template<typename inoutOperatorViewType,
890  typename fieldSignViewType>
892  inoutOperatorViewType _inoutOperator;
893  const fieldSignViewType _fieldSigns;
894 
895  KOKKOS_INLINE_FUNCTION
896  F_applyLeftFieldSigns( inoutOperatorViewType inoutOperator_,
897  fieldSignViewType fieldSigns_ )
898  : _inoutOperator(inoutOperator_), _fieldSigns(fieldSigns_) {}
899 
900  KOKKOS_INLINE_FUNCTION
901  void operator()(const ordinal_type cl) const {
902  const ordinal_type nlbf = _inoutOperator.extent(1);
903  const ordinal_type nrbf = _inoutOperator.extent(2);
904 
905  for (ordinal_type lbf=0;lbf<nlbf;++lbf)
906  for (ordinal_type rbf=0;rbf<nrbf;++rbf)
907  _inoutOperator(cl, lbf, rbf) *= _fieldSigns(cl, lbf);
908  }
909  };
910  }
911 
912  template<typename DeviceType>
913  template<typename inoutOperatorValueType, class ...inoutOperatorProperties,
914  typename fieldSignValueType, class ...fieldSignProperties>
915  void
917  applyLeftFieldSigns( Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...> inoutOperator,
918  const Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSigns ) {
919 
920 #ifdef HAVE_INTREPID2_DEBUG
921  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.rank() != 3, std::invalid_argument,
922  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Input operator container must have rank 3.");
923  INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
924  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Input field signs container must have rank 2.");
925  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(0) != fieldSigns.extent(0), std::invalid_argument,
926  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Zeroth dimensions (number of cells) of the operator and field signs containers must agree!");
927  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(1) != fieldSigns.extent(1), std::invalid_argument,
928  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): First dimensions (number of left fields) of the operator and field signs containers must agree!");
929 #endif
930 
931  constexpr bool are_accessible =
932  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
933  typename decltype(inoutOperator)::memory_space>::accessible &&
934  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
935  typename decltype(fieldSigns)::memory_space>::accessible;
936  static_assert(are_accessible, "FunctionSpaceTools<DeviceType>::applyLeftFieldSigns(..): input/output views' memory spaces are not compatible with DeviceType");
937 
939  <decltype(inoutOperator),decltype(fieldSigns)>;
940 
941  const ordinal_type C = inoutOperator.extent(0);
942  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
943  Kokkos::parallel_for( policy, FunctorType(inoutOperator, fieldSigns) );
944  }
945 
946  // ------------------------------------------------------------------------------------
947 
948  namespace FunctorFunctionSpaceTools {
952  template<typename inoutOperatorViewType,
953  typename fieldSignViewType>
955  inoutOperatorViewType _inoutOperator;
956  const fieldSignViewType _fieldSigns;
957 
958  KOKKOS_INLINE_FUNCTION
959  F_applyRightFieldSigns( inoutOperatorViewType inoutOperator_,
960  fieldSignViewType fieldSigns_ )
961  : _inoutOperator(inoutOperator_), _fieldSigns(fieldSigns_) {}
962 
963  KOKKOS_INLINE_FUNCTION
964  void operator()(const ordinal_type cl) const {
965  const ordinal_type nlbf = _inoutOperator.extent(1);
966  const ordinal_type nrbf = _inoutOperator.extent(2);
967 
968  for (ordinal_type lbf=0;lbf<nlbf;++lbf)
969  for (ordinal_type rbf=0;rbf<nrbf;++rbf)
970  _inoutOperator(cl, lbf, rbf) *= _fieldSigns(cl, rbf);
971  }
972  };
973  }
974 
975  template<typename DeviceType>
976  template<typename inoutOperatorValueType, class ...inoutOperatorProperties,
977  typename fieldSignValueType, class ...fieldSignProperties>
978  void
980  applyRightFieldSigns( Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...> inoutOperator,
981  const Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSigns ) {
982 
983 #ifdef HAVE_INTREPID2_DEBUG
984  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.rank() != 3, std::invalid_argument,
985  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Input operator container must have rank 3.");
986  INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
987  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Input field signs container must have rank 2.");
988  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(0) != fieldSigns.extent(0), std::invalid_argument,
989  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Zeroth dimensions (number of cells) of the operator and field signs containers must agree!");
990  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(2) != fieldSigns.extent(1), std::invalid_argument,
991  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Second dimension of the operator container and first dimension of the field signs container (number of right fields) must agree!");
992 #endif
993 
994  constexpr bool are_accessible =
995  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
996  typename decltype(inoutOperator)::memory_space>::accessible &&
997  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
998  typename decltype(fieldSigns)::memory_space>::accessible;
999  static_assert(are_accessible, "FunctionSpaceTools<DeviceType>::applyRightFieldSigns(..): input/output views' memory spaces are not compatible with DeviceType");
1000 
1002  <decltype(inoutOperator),decltype(fieldSigns)>;
1003 
1004  const ordinal_type C = inoutOperator.extent(0);
1005  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
1006  Kokkos::parallel_for( policy, FunctorType(inoutOperator, fieldSigns) );
1007  }
1008 
1009  // ------------------------------------------------------------------------------------
1010 
1011  namespace FunctorFunctionSpaceTools {
1015  template<typename inoutFunctionViewType,
1016  typename fieldSignViewType>
1018  inoutFunctionViewType _inoutFunction;
1019  const fieldSignViewType _fieldSigns;
1020 
1021  KOKKOS_INLINE_FUNCTION
1022  F_applyFieldSigns( inoutFunctionViewType inoutFunction_,
1023  fieldSignViewType fieldSigns_)
1024  : _inoutFunction(inoutFunction_), _fieldSigns(fieldSigns_) {}
1025 
1026  KOKKOS_INLINE_FUNCTION
1027  void operator()(const ordinal_type cl) const {
1028  const ordinal_type nbfs = _inoutFunction.extent(1);
1029  const ordinal_type npts = _inoutFunction.extent(2);
1030  const ordinal_type iend = _inoutFunction.extent(3);
1031  const ordinal_type jend = _inoutFunction.extent(4);
1032 
1033  for (ordinal_type bf=0;bf<nbfs;++bf)
1034  for (ordinal_type pt=0;pt<npts;++pt)
1035  for (ordinal_type i=0;i<iend;++i)
1036  for (ordinal_type j=0;j<jend;++j)
1037  _inoutFunction(cl, bf, pt, i, j) *= _fieldSigns(cl, bf);
1038  }
1039  };
1040  }
1041 
1042  template<typename DeviceType>
1043  template<typename inoutFunctionValueType, class ...inoutFunctionProperties,
1044  typename fieldSignValueType, class ...fieldSignProperties>
1045  void
1047  applyFieldSigns( Kokkos::DynRankView<inoutFunctionValueType,inoutFunctionProperties...> inoutFunction,
1048  const Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSigns ) {
1049 
1050 #ifdef HAVE_INTREPID2_DEBUG
1051  INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.rank() < 2 || inoutFunction.rank() > 5, std::invalid_argument,
1052  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Input function container must have rank 2, 3, 4, or 5.");
1053  INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
1054  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Input field signs container must have rank 2.");
1055  INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.extent(0) != fieldSigns.extent(0), std::invalid_argument,
1056  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Zeroth dimensions (number of integration domains) of the function and field signs containers must agree!");
1057  INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.extent(1) != fieldSigns.extent(1), std::invalid_argument,
1058  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): First dimensions (number of fields) of the function and field signs containers must agree!");
1059 
1060 #endif
1061 
1062  constexpr bool are_accessible =
1063  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1064  typename decltype(inoutFunction)::memory_space>::accessible &&
1065  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1066  typename decltype(fieldSigns)::memory_space>::accessible;
1067  static_assert(are_accessible, "FunctionSpaceTools<DeviceType>::applyFieldSigns(..): input/output views' memory spaces are not compatible with DeviceType");
1068 
1070  <decltype(inoutFunction),decltype(fieldSigns)>;
1071 
1072  const ordinal_type C = inoutFunction.extent(0);
1073  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
1074  Kokkos::parallel_for( policy, FunctorType(inoutFunction, fieldSigns) );
1075  }
1076 
1077  // ------------------------------------------------------------------------------------
1078 
1079  namespace FunctorFunctionSpaceTools {
1080 
1084  template<typename outputPointViewType,
1085  typename inputCoeffViewType,
1086  typename inputFieldViewType>
1088  outputPointViewType _outputPointVals;
1089  const inputCoeffViewType _inputCoeffs;
1090  const inputFieldViewType _inputFields;
1091 
1092  KOKKOS_INLINE_FUNCTION
1093  F_evaluateScalar( outputPointViewType outputPointVals_,
1094  inputCoeffViewType inputCoeffs_,
1095  inputFieldViewType inputFields_ )
1096  : _outputPointVals(outputPointVals_), _inputCoeffs(inputCoeffs_), _inputFields(inputFields_) {}
1097 
1098  KOKKOS_INLINE_FUNCTION
1099  void operator()(const ordinal_type cl, const ordinal_type pt) const {
1100  const ordinal_type nbfs = _inputFields.extent(1);
1101 
1102  for (ordinal_type bf=0;bf<nbfs;++bf)
1103  _outputPointVals(cl, pt) += _inputCoeffs(cl, bf) * _inputFields(cl, bf, pt);
1104  }
1105  };
1106 
1107  template<typename outputPointViewType,
1108  typename inputCoeffViewType,
1109  typename inputFieldViewType>
1111  outputPointViewType _outputPointVals;
1112  const inputCoeffViewType _inputCoeffs;
1113  const inputFieldViewType _inputFields;
1114 
1115  KOKKOS_INLINE_FUNCTION
1116  F_evaluateVector( outputPointViewType outputPointVals_,
1117  inputCoeffViewType inputCoeffs_,
1118  inputFieldViewType inputFields_ )
1119  : _outputPointVals(outputPointVals_), _inputCoeffs(inputCoeffs_), _inputFields(inputFields_) {}
1120 
1121  KOKKOS_INLINE_FUNCTION
1122  void operator()(const ordinal_type cl, const ordinal_type pt) const {
1123  const ordinal_type nbfs = _inputFields.extent(1);
1124  const ordinal_type iend = _inputFields.extent(3);
1125 
1126  for (ordinal_type bf=0;bf<nbfs;++bf)
1127  for (ordinal_type i=0;i<iend;++i)
1128  _outputPointVals(cl, pt, i) += _inputCoeffs(cl, bf) * _inputFields(cl, bf, pt, i);
1129  }
1130  };
1131 
1132  template<typename outputPointViewType,
1133  typename inputCoeffViewType,
1134  typename inputFieldViewType>
1136  outputPointViewType _outputPointVals;
1137  const inputCoeffViewType _inputCoeffs;
1138  const inputFieldViewType _inputFields;
1139 
1140  KOKKOS_INLINE_FUNCTION
1141  F_evaluateTensor( outputPointViewType outputPointVals_,
1142  inputCoeffViewType inputCoeffs_,
1143  inputFieldViewType inputFields_ )
1144  : _outputPointVals(outputPointVals_), _inputCoeffs(inputCoeffs_), _inputFields(inputFields_) {}
1145 
1146  KOKKOS_INLINE_FUNCTION
1147  void operator()(const ordinal_type cl, const ordinal_type pt) const {
1148  const ordinal_type nbfs = _inputFields.extent(1);
1149  const ordinal_type iend = _inputFields.extent(3);
1150  const ordinal_type jend = _inputFields.extent(4);
1151 
1152  for (ordinal_type bf=0;bf<nbfs;++bf)
1153  for (ordinal_type i=0;i<iend;++i)
1154  for (ordinal_type j=0;j<jend;++j)
1155  _outputPointVals(cl, pt, i, j) += _inputCoeffs(cl, bf) * _inputFields(cl, bf, pt, i, j);
1156  }
1157  };
1158  }
1159 
1160  template<typename DeviceType>
1161  template<typename outputPointValueType, class ...outputPointProperties,
1162  typename inputCoeffValueType, class ...inputCoeffProperties,
1163  typename inputFieldValueType, class ...inputFieldProperties>
1164  void
1166  evaluate( Kokkos::DynRankView<outputPointValueType,outputPointProperties...> outputPointVals,
1167  const Kokkos::DynRankView<inputCoeffValueType, inputCoeffProperties...> inputCoeffs,
1168  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields ) {
1169 
1170 #ifdef HAVE_INTREPID2_DEBUG
1171  INTREPID2_TEST_FOR_EXCEPTION( inputFields.rank() < 3 || inputFields.rank() > 5, std::invalid_argument,
1172  ">>> ERROR (FunctionSpaceTools::evaluate): Input fields container must have rank 3, 4, or 5.");
1173  INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.rank() != 2, std::invalid_argument,
1174  ">>> ERROR (FunctionSpaceTools::evaluate): Input coefficient container must have rank 2.");
1175  INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.rank() != (inputFields.rank()-1), std::invalid_argument,
1176  ">>> ERROR (FunctionSpaceTools::evaluate): Output values container must have rank one less than the rank of the input fields container.");
1177  INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.extent(0) != inputFields.extent(0), std::invalid_argument,
1178  ">>> ERROR (FunctionSpaceTools::evaluate): Zeroth dimensions (number of cells) of the coefficient and fields input containers must agree!");
1179  INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.extent(1) != inputFields.extent(1), std::invalid_argument,
1180  ">>> ERROR (FunctionSpaceTools::evaluate): First dimensions (number of fields) of the coefficient and fields input containers must agree!");
1181  INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.extent(0) != inputFields.extent(0), std::invalid_argument,
1182  ">>> ERROR (FunctionSpaceTools::evaluate): Zeroth dimensions (number of cells) of the input fields container and the output values container must agree!");
1183  for (size_type i=1;i<outputPointVals.rank();++i)
1184  INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.extent(i) != inputFields.extent(i+1), std::invalid_argument,
1185  ">>> ERROR (FunctionSpaceTools::evaluate): outputPointVals dimension(i) does not match to inputFields dimension(i+1).");
1186 #endif
1187 
1188  constexpr bool are_accessible =
1189  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1190  typename decltype(outputPointVals)::memory_space>::accessible &&
1191  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1192  typename decltype(inputCoeffs)::memory_space>::accessible &&
1193  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1194  typename decltype(inputFields)::memory_space>::accessible;
1195  static_assert(are_accessible, "FunctionSpaceTools<DeviceType>::evaluate(..): input/output views' memory spaces are not compatible with DeviceType");
1196 
1197  using range_policy_type = Kokkos::MDRangePolicy< ExecSpaceType, Kokkos::Rank<2>, Kokkos::IndexType<ordinal_type> >;
1198 
1199  const range_policy_type policy( { 0, 0 },
1200  { /*C*/ inputFields.extent(0), /*P*/ inputFields.extent(2)} );
1201 
1202  if (inputFields.rank() == 3) {
1204  Kokkos::parallel_for( policy, FunctorType(outputPointVals, inputCoeffs, inputFields) );
1205  }
1206  else if (inputFields.rank() == 4) {
1208  Kokkos::parallel_for( policy, FunctorType(outputPointVals, inputCoeffs, inputFields) );
1209  }
1210  else {
1212  Kokkos::parallel_for( policy, FunctorType(outputPointVals, inputCoeffs, inputFields) );
1213  }
1214  }
1215 
1216  // ------------------------------------------------------------------------------------
1217 
1218 } // end namespace Intrepid2
1219 
1220 #endif
static void vectorNorm(Kokkos::DynRankView< normArrayValueType, normArrayProperties...> normArray, const Kokkos::DynRankView< inVecValueType, inVecProperties...> inVecs, const ENorm normType)
Computes norms (1, 2, infinity) of vectors stored in a array of total rank 2 (array of vectors)...
static void mapHDivDataFromPhysToRef(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianInverseValueType, jacobianInverseProperties...> jacobianInv, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of a (vector) data in the H-div space, defined in the physical space, stored in the user-provided container inputVals and indexed by (C,P,D), into the output container outputVals, defined on the reference cell and indexed by (C,P,D).
Functor for applyLeftFieldSigns, see Intrepid2::FunctionSpaceTools for more.
static void contractFieldFieldTensor(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< leftFieldValueType, leftFieldProperties...> leftFields, const Kokkos::DynRankView< rightFieldValueType, rightFieldProperties...> rightFields, const bool sumInto=false)
Contracts the &quot;point&quot; and &quot;space&quot; dimensions P, D1, and D2 of two rank-5 containers with dimensions (...
static void mapHGradDataFromPhysSideToRefSide(Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input)
Transformation of a (scalar) data in the H-grad space, defined in physical space, stored in the user-...
static void applyFieldSigns(Kokkos::DynRankView< inoutFunctionValueType, inoutFunctionProperties...> inoutFunction, const Kokkos::DynRankView< fieldSignValueType, fieldSignProperties...> fieldSigns)
Applies field signs, stored in the user-provided container fieldSigns and indexed by (C...
static void scalarMultiplyDataData(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight, const bool reciprocal=false)
There are two use cases: (1) multiplies a rank-2, 3, or 4 container inputDataRight with dimensions (C...
static void contractFieldFieldScalar(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< leftFieldValueType, leftFieldProperties...> leftFields, const Kokkos::DynRankView< rightFieldValueType, rightFieldProperties...> rightFields, const bool sumInto=false)
Contracts the &quot;point&quot; dimension P of two rank-3 containers with dimensions (C,L,P) and (C...
static void HVOLtransformVALUE(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of a (scalar) value field in the H-vol space, defined at points on a reference cell...
static void matmatProductDataData(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight, const char transpose= 'N')
There are two use cases: (1) matrix-matrix product of a rank-4 container inputDataRight with dimensio...
static void applyRightFieldSigns(Kokkos::DynRankView< inoutOperatorValueType, inoutOperatorProperties...> inoutOperator, const Kokkos::DynRankView< fieldSignValueType, fieldSignProperties...> fieldSigns)
Applies right (column) signs, stored in the user-provided container fieldSigns and indexed by (C...
static void vectorMultiplyDataField(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields)
Cross or outer product of data and fields; please read the description below.
static void contractDataFieldTensor(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields, const bool sumInto=false)
Contracts the &quot;point&quot; and &quot;space&quot; dimensions P, D1 and D2 of a rank-5 container and a rank-4 containe...
static void dotMultiplyDataField(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields)
Dot product of data and fields; please read the description below.
static void crossProductDataData(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight)
There are two use cases: (1) cross product of a rank-3 container inputDataRight with dimensions (C...
Defines TensorArgumentIterator, which allows systematic enumeration of a TensorData object...
static void crossProductDataField(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields)
There are two use cases: (1) cross product of a rank-4 container inputFields with dimensions (C...
static void clone(Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input)
Clone input array.
Functor for applyRightFieldSigns, see Intrepid2::FunctionSpaceTools for more.
static void contractDataDataScalar(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight, const bool sumInto=false)
Contracts the &quot;point&quot; dimensions P of rank-2 containers with dimensions (C,P), and returns the result...
static void evaluate(Kokkos::DynRankView< outputPointValueType, outputPointProperties...> outputPointVals, const Kokkos::DynRankView< inputCoeffValueType, inputCoeffProperties...> inputCoeffs, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields)
Computes point values outPointVals of a discrete function specified by the basis inFields and coeffic...
static void matvecProductDataField(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields, const char transpose= 'N')
There are two use cases: (1) matrix-vector product of a rank-4 container inputFields with dimensions ...
static void scalarMultiplyDataData(Kokkos::DynRankView< outputDataValuetype, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight, const bool reciprocal=false)
Scalar multiplication of data and data; please read the description below.
static void dotMultiplyDataData(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight)
There are two use cases: (1) dot product of a rank-2, 3 or 4 container inputDataRight with dimensions...
static void outerProductDataField(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields)
There are two use cases: (1) outer product of a rank-4 container inputFields with dimensions (C...
static void cloneFields(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields)
Replicates a rank-2, 3, or 4 container with dimensions (F,P), (F,P,D1) or (F,P,D1,D2), representing the values of a scalar, vector or a tensor field, into an output value container of size (C,F,P), (C,F,P,D1) or (C,F,P,D1,D2).
static void scalarMultiplyDataField(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields, const bool reciprocal=false)
There are two use cases: (1) multiplies a rank-3, 4, or 5 container inputFields with dimensions (C...
static void contractDataFieldVector(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields, const bool sumInto=false)
Contracts the &quot;point&quot; and &quot;space&quot; dimensions P and D of a rank-4 container and a rank-3 container wit...
static void contractDataFieldScalar(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields, const bool sumInto=false)
Contracts the &quot;point&quot; dimensions P of a rank-3 containers and a rank-2 container with dimensions (C...
static void mapHCurlDataCrossNormalFromPhysSideToRefSide(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< tangentsValueType, tangentsProperties...> tangents, const Kokkos::DynRankView< metricTensorInvValueType, metricTensorInvProperties...> metricTensorInv, const Kokkos::DynRankView< metricTensorDetValueType, metricTensorDetProperties...> metricTensorDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of 3D HCURL data from physical side to reference side. It takes the input vector defi...
static void dotMultiplyDataData(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight)
Dot product of data and data; please read the description below.
static void contractDataDataVector(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight, const bool sumInto=false)
Contracts the &quot;point&quot; and &quot;space&quot; dimensions P and D of rank-3 containers with dimensions (C...
static void matmatProductDataField(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields, const char transpose= 'N')
There are two use cases: (1) matrix-matrix product of a rank-5 container inputFields with dimensions ...
static void vectorMultiplyDataData(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight)
Cross or outer product of data and data; please read the description below.
static void multiplyMeasure(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< inputMeasureValueType, inputMeasureProperties...> inputMeasure, const Kokkos::DynRankView< inputValValueType, inputValProperteis...> inputVals)
Multiplies fields inputVals by weighted measures inputMeasure and returns the field array outputVals;...
static void mapHGradDataFromPhysToRef(Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input)
Transformation of a (scalar) data in the H-grad space, defined in physical space, stored in the user-...
Defines the Intrepid2::FunctorIterator class, as well as the Intrepid2::functor_returns_ref SFINAE he...
static void outerProductDataData(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValuetype, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight)
There are two use cases: (1) outer product of a rank-3 container inputDataRight with dimensions (C...
static void computeEdgeMeasure(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< inputJacValueType, inputJacProperties...> inputJac, const Kokkos::DynRankView< inputWeightValueType, inputWeightProperties...> inputWeights, const int whichEdge, const shards::CellTopology parentCell, const Kokkos::DynRankView< scratchValueType, scratchProperties...> scratch)
Returns the weighted integration measures outVals with dimensions (C,P) used for the computation of e...
Functor to evaluate functions, see Intrepid2::FunctionSpaceTools for more.
static void integrate(Kokkos::DynRankView< outputValueValueType, outputValueProperties...> outputValues, const Kokkos::DynRankView< leftValueValueType, leftValueProperties...> leftValues, const Kokkos::DynRankView< rightValueValueType, rightValueProperties...> rightValues, const bool sumInto=false)
Contracts leftValues and rightValues arrays on the point and possibly space dimensions and stores the...
static void mapHVolDataFromPhysToRef(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of a (scalar) data in the H-vol space, defined in the physical space, stored in the user-provided container inputVals and indexed by (C,P), into the output container outputVals, defined on the reference cell and indexed by (C,P).
static void tensorMultiplyDataData(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight, const char transpose= 'N')
Matrix-vector or matrix-matrix product of data and data; please read the description below...
static void HGRADtransformGRAD(OutputValViewType outputVals, const JacobianInverseViewType jacobianInverse, const InputValViewType inputVals)
Transformation of a gradient field in the H-grad space, defined at points on a reference cell...
static void computeFaceMeasure(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< inputJacValueType, inputJacProperties...> inputJac, const Kokkos::DynRankView< inputWeightValueType, inputWeightPropertes...> inputWeights, const int whichFace, const shards::CellTopology parentCell, const Kokkos::DynRankView< scratchValueType, scratchProperties...> scratch)
Returns the weighted integration measures outputVals with dimensions (C,P) used for the computation o...
static void applyLeftFieldSigns(Kokkos::DynRankView< inoutOperatorValueType, inoutOperatorProperties...> inoutOperator, const Kokkos::DynRankView< fieldSignValueType, fieldSignProperties...> fieldSigns)
Applies left (row) signs, stored in the user-provided container fieldSigns and indexed by (C...
static void mapHDivDataDotNormalFromPhysSideToRefSide(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> metricTensorDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of HDIV data from physical side to reference side. It takes the input defined on phys...
static void tensorMultiplyDataField(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields, const char transpose= 'N')
Matrix-vector or matrix-matrix product of data and fields; please read the description below...
static void contractDataDataTensor(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight, const bool sumInto=false)
Contracts the &quot;point&quot; and &quot;space&quot; dimensions P, D1 and D2 of rank-4 containers with dimensions (C...
static void scalarMultiplyDataField(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataPropertes...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields, const bool reciprocal=false)
Scalar multiplication of data and fields; please read the description below.
static void HCURLtransformCURL(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianValueType, jacobianProperties...> jacobian, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of a 3D curl field in the H-curl space, defined at points on a reference cell...
static void getPhysicalFaceNormals(Kokkos::DynRankView< faceNormalValueType, faceNormalProperties...> faceNormals, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetFaceOrd, const shards::CellTopology parentCell)
Computes non-normalized normal vectors to physical faces in a face workset ; (see Subcell worksets fo...
static void HCURLtransformVALUE(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianInverseValueType, jacobianInverseProperties...> jacobianInverse, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of a (vector) value field in the H-curl space, defined at points on a reference cell...
static void HDIVtransformVALUE(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianValueType, jacobianProperties...> jacobian, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of a (vector) value field in the H-div space, defined at points on a reference cell...
static void HDIVtransformDIV(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of a divergence field in the H-div space, defined at points on a reference cell...
Functor for calculation of cell measure, see Intrepid2::FunctionSpaceTools for more.
static void HGRADtransformCURL(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianValueType, jacobianProperties...> jacobian, const Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties...> jacobianDet, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of a 2D curl field in the H-grad space, defined at points on a reference cell...
static void contractFieldFieldVector(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< leftFieldValueType, leftFieldProperties...> leftFields, const Kokkos::DynRankView< rightFieldValueType, rightFieldProperties...> rightFields, const bool sumInto=false)
Contracts the &quot;point&quot; and &quot;space&quot; dimensions P and D1 of two rank-4 containers with dimensions (C...
static bool computeCellMeasure(OutputValViewType outputVals, const InputDetViewType inputDet, const InputWeightViewType inputWeights)
Returns the weighted integration measures outputVals with dimensions (C,P) used for the computation o...
static void getPhysicalEdgeTangents(Kokkos::DynRankView< edgeTangentValueType, edgeTangentProperties...> edgeTangents, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetEdgeOrd, const shards::CellTopology parentCell)
Computes non-normalized tangent vectors to physical edges in an edge workset ; (see Subcell worksets ...
static void dotMultiplyDataField(Kokkos::DynRankView< outputFieldValueType, outputFieldProperties...> outputFields, const Kokkos::DynRankView< inputDataValueType, inputDataProperties...> inputData, const Kokkos::DynRankView< inputFieldValueType, inputFieldProperties...> inputFields)
There are two use cases: (1) dot product of a rank-3, 4 or 5 container inputFields with dimensions (C...
Functor for applyFieldSigns, see Intrepid2::FunctionSpaceTools for more.
static void HGRADtransformVALUE(Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input)
Transformation of a (scalar) value field in the H-grad space, defined at points on a reference cell...
static void mapHCurlDataFromPhysToRef(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianValueType, jacobianProperties...> jacobian, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of a (vector) data in the H-curl space, defined in the physical space, stored in the user-provided container inputVals and indexed by (C,P,D), into the output container outputVals, defined on the reference cell and indexed by (C,P,D).
static void matvecProductDataData(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight, const char transpose= 'N')
There are two use cases: (1) matrix-vector product of a rank-3 container inputDataRight with dimensio...