Intrepid2
Intrepid2_FunctionSpaceToolsDef.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid2 Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Kyungjoo Kim (kyukim@sandia.gov), or
38 // Mauro Perego (mperego@sandia.gov)
39 //
40 // ************************************************************************
41 // @HEADER
42 
49 #ifndef __INTREPID2_FUNCTIONSPACETOOLS_DEF_HPP__
50 #define __INTREPID2_FUNCTIONSPACETOOLS_DEF_HPP__
51 
52 namespace Intrepid2 {
53 
54  // ------------------------------------------------------------------------------------
55  template<typename SpT>
56  template<typename outputValueType, class ...outputProperties,
57  typename inputValueType, class ...inputProperties>
58  void
60  HGRADtransformVALUE( Kokkos::DynRankView<outputValueType,outputProperties...> output,
61  const Kokkos::DynRankView<inputValueType, inputProperties...> input ) {
62  if(output.rank() == input.rank()) {
63 #ifdef HAVE_INTREPID2_DEBUG
64  {
65  for (size_type i=0;i< input.rank();++i) {
66  INTREPID2_TEST_FOR_EXCEPTION( (input.extent(i) != output.extent(i)), std::invalid_argument,
67  ">>> ERROR (FunctionSpaceTools::HGRADtransformVALUE): Dimensions of input and output fields containers do not match.");
68  }
69  }
70 #endif
71  RealSpaceTools<SpT>::clone(output, input);
72  }
73  else
74  ArrayTools<SpT>::cloneFields(output, input);
75  }
76 
77  // ------------------------------------------------------------------------------------
78 
79  namespace FunctorFunctionSpaceTools {
83  template <typename OutputViewType,
84  typename jacInverseViewType,
85  typename inputViewType,
86  ordinal_type spaceDim>
88  OutputViewType _output;
89  const jacInverseViewType _jacInverse;
90  const inputViewType _input;
91 
92  // output CPDD, left CPDD or PDD, right FPD
93  KOKKOS_INLINE_FUNCTION
94  F_HGRADtransformGRAD(OutputViewType output_,
95  jacInverseViewType jacInverse_,
96  inputViewType input_)
97  : _output(output_),
98  _jacInverse(jacInverse_),
99  _input(input_) {}
100 
101  KOKKOS_INLINE_FUNCTION
102  void operator()(const ordinal_type cl,
103  const ordinal_type bf,
104  const ordinal_type pt) const {
105  /* */ auto y = Kokkos::subview(_output, cl, bf, pt, Kokkos::ALL());
106  const auto A = Kokkos::subview(_jacInverse, cl, pt, Kokkos::ALL(), Kokkos::ALL());
107  const auto x = Kokkos::subview(_input, bf, pt, Kokkos::ALL());
108 
109  if (spaceDim == 2) {
110  Kernels::Serial::matvec_trans_product_d2( y, A, x );
111  } else {
112  Kernels::Serial::matvec_trans_product_d3( y, A, x );
113  }
114  }
115  };
116  }
117 
118  template<typename SpT>
119  template<typename outputValValueType, class ...outputValProperties,
120  typename jacobianInverseValueType, class ...jacobianInverseProperties,
121  typename inputValValueType, class ...inputValProperties>
122  void
124  HGRADtransformGRAD( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
125  const Kokkos::DynRankView<jacobianInverseValueType,jacobianInverseProperties...> jacobianInverse,
126  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
127  return HCURLtransformVALUE(outputVals, jacobianInverse, inputVals);
128 
129  // this modification is for 2d and 3d (not 1d)
130  // this is an attempt to measure the overhead of subview of dynrankview.
131 
132  // typedef Kokkos::DynRankView<outputValValueType,outputValProperties...> OutputViewType;
133  // typedef const Kokkos::DynRankView<jacobianInverseValueType,jacobianInverseProperties...> jacInverseViewType;
134  // typedef const Kokkos::DynRankView<inputValValueType,inputValProperties...> inputViewType;
135 
136  // const ordinal_type
137  // C = outputVals.extent(0),
138  // F = outputVals.extent(1),
139  // P = outputVals.extent(2);
140 
141  // using range_policy_type = Kokkos::Experimental::MDRangePolicy
142  // < SpT, Kokkos::Experimental::Rank<3>, Kokkos::IndexType<ordinal_type> >;
143  // range_policy_type policy( { 0, 0, 0 },
144  // { C, F, P } );
145 
146  // const ordinal_type spaceDim = inputVals.extent(2);
147  // switch (spaceDim) {
148  // case 2: {
149  // typedef FunctorFunctionSpaceTools::F_HGRADtransformGRAD<OutputViewType, jacInverseViewType, inputViewType, 2> FunctorType;
150  // Kokkos::parallel_for( policy, FunctorType(outputVals, jacobianInverse, inputVals) );
151  // break;
152  // }
153  // case 3: {
154  // typedef FunctorFunctionSpaceTools::F_HGRADtransformGRAD<OutputViewType, jacInverseViewType, inputViewType, 3> FunctorType;
155  // Kokkos::parallel_for( policy, FunctorType(outputVals, jacobianInverse, inputVals) );
156  // break;
157  // }
158  // default: {
159  // INTREPID2_TEST_FOR_EXCEPTION( true, std::invalid_argument,
160  // ">>> ERROR (FunctionSpaceTools::HGRADtransformGRAD): spaceDim is not 2 or 3.");
161  // break;
162  // }
163  // }
164  }
165 
166  // ------------------------------------------------------------------------------------
167 
168  template<typename SpT>
169  template<typename outputValValueType, class ...outputValProperties,
170  typename jacobianInverseValueType, class ...jacobianInverseProperties,
171  typename inputValValueType, class ...inputValProperties>
172  void
174  HCURLtransformVALUE( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
175  const Kokkos::DynRankView<jacobianInverseValueType,jacobianInverseProperties...> jacobianInverse,
176  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
177  ArrayTools<SpT>::matvecProductDataField(outputVals, jacobianInverse, inputVals, 'T');
178  }
179 
180  // ------------------------------------------------------------------------------------
181 
182  template<typename SpT>
183  template<typename outputValValueType, class ...outputValProperties,
184  typename jacobianValueType, class ...jacobianProperties,
185  typename jacobianDetValueType, class ...jacobianDetProperties,
186  typename inputValValueType, class ...inputValProperties>
187  void
189  HCURLtransformCURL( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
190  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
191  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
192  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
193  if(jacobian.data()==NULL || jacobian.extent(2)==2) //2D case
194  return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
195  else
196  return HDIVtransformVALUE(outputVals, jacobian, jacobianDet, inputVals);
197  }
198 
199  // ------------------------------------------------------------------------------------
200 
201  template<typename SpT>
202  template<typename outputValValueType, class ...outputValProperties,
203  typename jacobianDetValueType, class ...jacobianDetProperties,
204  typename inputValValueType, class ...inputValProperties>
205  void
207  HCURLtransformCURL( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
208  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
209  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
210 #ifdef HAVE_INTREPID2_DEBUG
211  {
212  INTREPID2_TEST_FOR_EXCEPTION( outputVals.rank() == 4, std::invalid_argument,
213  ">>> ERROR (FunctionSpaceTools::HCURLtransformCURL): Output rank must have rank 3.\n If these are 3D fields, then use the appropriate overload of this function.");
214  }
215 #endif
216  return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
217  }
218 
219  // ------------------------------------------------------------------------------------
220 
221  template<typename SpT>
222  template<typename outputValValueType, class ...outputValProperties,
223  typename jacobianValueType, class ...jacobianProperties,
224  typename jacobianDetValueType, class ...jacobianDetProperties,
225  typename inputValValueType, class ...inputValProperties>
226  void
228  HGRADtransformCURL( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
229  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
230  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
231  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
232 #ifdef HAVE_INTREPID2_DEBUG
233  {
234  INTREPID2_TEST_FOR_EXCEPTION( outputVals.extent(3)!=2, std::invalid_argument,
235  ">>> ERROR (FunctionSpaceTools::HGRADtransformCURL):\n output field is 3D by the function is meant for 2D fields");
236  }
237 #endif
238  return HDIVtransformVALUE(outputVals, jacobian, jacobianDet, inputVals);
239  }
240 
241  // ------------------------------------------------------------------------------------
242 
243  template<typename SpT>
244  template<typename outputValValueType, class ...outputValProperties,
245  typename jacobianValueType, class ...jacobianProperties,
246  typename jacobianDetValueType, class ...jacobianDetProperties,
247  typename inputValValueType, class ...inputValProperties>
248  void
250  HDIVtransformVALUE( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
251  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
252  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
253  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
254  ArrayTools<SpT>::matvecProductDataField(outputVals, jacobian, inputVals, 'N');
255  ArrayTools<SpT>::scalarMultiplyDataField(outputVals, jacobianDet, outputVals, true);
256  }
257 
258  // ------------------------------------------------------------------------------------
259 
260  template<typename SpT>
261  template<typename outputValValueType, class ...outputValProperties,
262  typename jacobianDetValueType, class ...jacobianDetProperties,
263  typename inputValValueType, class ...inputValProperties>
264  void
266  HDIVtransformDIV( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
267  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
268  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
269  return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
270  }
271 
272  // ------------------------------------------------------------------------------------
273 
274  template<typename SpT>
275  template<typename outputValValueType, class ...outputValProperties,
276  typename jacobianDetValueType, class ...jacobianDetProperties,
277  typename inputValValueType, class ...inputValProperties>
278  void
280  HVOLtransformVALUE( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
281  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
282  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
283  ArrayTools<SpT>::scalarMultiplyDataField(outputVals, jacobianDet, inputVals, true);
284  }
285 
286  // ------------------------------------------------------------------------------------
287 
288  template<typename SpT>
289  template<typename outputValueValueType, class ...outputValueProperties,
290  typename leftValueValueType, class ...leftValueProperties,
291  typename rightValueValueType, class ...rightValueProperties>
292  void
294  integrate( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
295  const Kokkos::DynRankView<leftValueValueType, leftValueProperties...> leftValues,
296  const Kokkos::DynRankView<rightValueValueType, rightValueProperties...> rightValues,
297  const bool sumInto ) {
298 
299 #ifdef HAVE_INTREPID2_DEBUG
300  {
301  INTREPID2_TEST_FOR_EXCEPTION( leftValues.rank() < 2 ||
302  leftValues.rank() > 4, std::invalid_argument,
303  ">>> ERROR (FunctionSpaceTools::integrate): Left data must have rank 2, 3 or 4.");
304  INTREPID2_TEST_FOR_EXCEPTION( outputValues.rank() < 1 ||
305  outputValues.rank() > 3, std::invalid_argument,
306  ">>> ERROR (FunctionSpaceTools::integrate): Output values must have rank 1, 2 or 3.");
307  }
308 #endif
309 
310  const ordinal_type outRank = outputValues.rank();
311  const ordinal_type leftRank = leftValues.rank();
312  const ordinal_type mode = outRank*10 + leftRank;
313 
314  switch (mode) {
315  // DataData
316  case 12:
318  leftValues,
319  rightValues,
320  sumInto );
321  break;
322  case 13:
324  leftValues,
325  rightValues,
326  sumInto );
327  break;
328  case 14:
330  leftValues,
331  rightValues,
332  sumInto );
333  break;
334 
335  // DataField
336  case 22:
338  leftValues,
339  rightValues,
340  sumInto );
341  break;
342  case 23:
344  leftValues,
345  rightValues,
346  sumInto );
347  break;
348  case 24:
350  leftValues,
351  rightValues,
352  sumInto );
353  break;
354 
355  // FieldField
356  case 33:
358  leftValues,
359  rightValues,
360  sumInto );
361  break;
362  case 34:
364  leftValues,
365  rightValues,
366  sumInto );
367  break;
368  case 35:
370  leftValues,
371  rightValues,
372  sumInto );
373  break;
374  default: {
375  INTREPID2_TEST_FOR_EXCEPTION( outRank < 1 || outRank > 3, std::runtime_error,
376  ">>> ERROR (FunctionSpaceTools::integrate): outRank must be 1,2, or 3.");
377  INTREPID2_TEST_FOR_EXCEPTION( leftRank < 2 || leftRank > 4, std::runtime_error,
378  ">>> ERROR (FunctionSpaceTools::integrate): leftRank must be 1,2, 3 or 4.");
379  }
380  }
381  }
382 
383  // ------------------------------------------------------------------------------------
384  namespace FunctorFunctionSpaceTools {
388  template<typename outputValViewType,
389  typename inputDetViewType,
390  typename inputWeightViewType>
392  outputValViewType _outputVals;
393  const inputDetViewType _inputDet;
394  const inputWeightViewType _inputWeight;
395 
396  KOKKOS_INLINE_FUNCTION
397  F_computeCellMeasure( outputValViewType outputVals_,
398  inputDetViewType inputDet_,
399  inputWeightViewType inputWeight_)
400  : _outputVals(outputVals_),
401  _inputDet(inputDet_),
402  _inputWeight(inputWeight_) {}
403 
404  typedef ordinal_type value_type;
405 
406 // KOKKOS_INLINE_FUNCTION
407 // void init( value_type &dst ) const {
408 // dst = false;
409 // }
410 
411 // KOKKOS_INLINE_FUNCTION
412 // void join( volatile value_type &dst,
413 // const volatile value_type &src ) const {
414 // dst |= src;
415 // }
416 
417  KOKKOS_INLINE_FUNCTION
418  void operator()(const size_type cl, value_type &dst) const {
419  // negative jacobian check
420  const bool hasNegativeDet = (_inputDet(cl, 0) < 0.0);
421  dst |= hasNegativeDet;
422 
423  // make it absolute
424  const auto sign = (hasNegativeDet ? -1.0 : 1.0);
425  const ordinal_type pt_end = _outputVals.extent(1);
426  for (ordinal_type pt=0;pt<pt_end;++pt)
427  _outputVals(cl, pt) = sign*_inputDet(cl, pt)*_inputWeight(pt);
428  }
429  };
430  }
431 
432  template<typename SpT>
433  template<typename outputValValueType, class ...outputValProperties,
434  typename inputDetValueType, class ...inputDetProperties,
435  typename inputWeightValueType, class ...inputWeightProperties>
436  bool
438  computeCellMeasure( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
439  const Kokkos::DynRankView<inputDetValueType, inputDetProperties...> inputDet,
440  const Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeights ) {
441 #ifdef HAVE_INTREPID2_DEBUG
442  {
443  INTREPID2_TEST_FOR_EXCEPTION( inputDet.rank() != 2 ||
444  inputWeights.rank() != 1 ||
445  outputVals.rank() != 2, std::invalid_argument,
446  ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Ranks are not compatible.");
447  INTREPID2_TEST_FOR_EXCEPTION( outputVals.extent(0) != inputDet.extent(0), std::invalid_argument,
448  ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Cell dimension does not match.");
449  INTREPID2_TEST_FOR_EXCEPTION( inputDet.extent(1) != outputVals.extent(1) ||
450  inputWeights.extent(0) != outputVals.extent(1), std::invalid_argument,
451  ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Point dimension does not match.");
452  }
453 #endif
454  typedef Kokkos::DynRankView<outputValValueType, outputValProperties...> outputValViewType;
455  typedef Kokkos::DynRankView<inputDetValueType, inputDetProperties...> inputDetViewType;
456  typedef Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeightViewType;
458  <outputValViewType,inputDetViewType,inputWeightViewType> FunctorType;
459 
460  const ordinal_type C = inputDet.extent(0);
461  Kokkos::RangePolicy<SpT,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
462 
463  typename FunctorType::value_type hasNegativeDet = false;
464  Kokkos::parallel_reduce( policy, FunctorType(outputVals, inputDet, inputWeights), hasNegativeDet );
465 
466  return hasNegativeDet;
467  }
468 
469  // ------------------------------------------------------------------------------------
470 
471  template<typename SpT>
472  template<typename outputValValueType, class ...outputValProperties,
473  typename inputJacValueType, class ...inputJacProperties,
474  typename inputWeightValueType, class ...inputWeightProperties,
475  typename scratchValueType, class ...scratchProperties>
476  void
478  computeFaceMeasure( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
479  const Kokkos::DynRankView<inputJacValueType, inputJacProperties...> inputJac,
480  const Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeights,
481  const ordinal_type whichFace,
482  const shards::CellTopology parentCell,
483  const Kokkos::DynRankView<scratchValueType, scratchProperties...> scratch ) {
484 #ifdef HAVE_INTREPID2_DEBUG
485  INTREPID2_TEST_FOR_EXCEPTION( inputJac.rank() != 4, std::invalid_argument,
486  ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Input Jacobian container must have rank 4.");
487  INTREPID2_TEST_FOR_EXCEPTION( scratch.rank() != 1, std::invalid_argument,
488  ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Scratch view imust have rank 1.");
489  INTREPID2_TEST_FOR_EXCEPTION( scratch.span() < inputJac.span(), std::invalid_argument,
490  ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Scratch storage must be greater than or equal to inputJac's one.");
491 #endif
492 
493  // face normals (reshape scratch)
494  // Kokkos::DynRankView<scratchValueType,scratchProperties...> faceNormals(scratch.data(),
495  // inputJac.extent(0),
496  // inputJac.extent(1),
497  // inputJac.extent(2));
498  auto vcprop = Kokkos::common_view_alloc_prop(scratch);
499  //typedef Kokkos::DynRankView<scratchValueType, typename decltype(scratch)::memory_space> viewType;
500  typedef Kokkos::DynRankView<scratchValueType, SpT> viewType;
501  viewType faceNormals(Kokkos::view_wrap(scratch.data(), vcprop),
502  inputJac.extent(0),
503  inputJac.extent(1),
504  inputJac.extent(2));
505 
506  // compute normals
507  CellTools<SpT>::getPhysicalFaceNormals(faceNormals, inputJac, whichFace, parentCell);
508 
509  // compute lenghts of normals
510  RealSpaceTools<SpT>::vectorNorm(outputVals, faceNormals, NORM_TWO);
511 
512  // multiply with weights
513  ArrayTools<SpT>::scalarMultiplyDataData(outputVals, outputVals, inputWeights);
514  }
515 
516  // ------------------------------------------------------------------------------------
517 
518  template<typename SpT>
519  template<typename outputValValueType, class ...outputValProperties,
520  typename inputJacValueType, class ...inputJacProperties,
521  typename inputWeightValueType, class ...inputWeightProperties,
522  typename scratchValueType, class ...scratchProperties>
523  void
525  computeEdgeMeasure( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
526  const Kokkos::DynRankView<inputJacValueType, inputJacProperties...> inputJac,
527  const Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeights,
528  const ordinal_type whichEdge,
529  const shards::CellTopology parentCell,
530  const Kokkos::DynRankView<scratchValueType, scratchProperties...> scratch ) {
531 #ifdef HAVE_INTREPID2_DEBUG
532  INTREPID2_TEST_FOR_EXCEPTION( (inputJac.rank() != 4), std::invalid_argument,
533  ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Input Jacobian container must have rank 4.");
534  INTREPID2_TEST_FOR_EXCEPTION( scratch.rank() != 1, std::invalid_argument,
535  ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Scratch view must have a rank 1.");
536  INTREPID2_TEST_FOR_EXCEPTION( scratch.span() < inputJac.span(), std::invalid_argument,
537  ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Scratch storage must be greater than or equal to inputJac'one.");
538 #endif
539 
540  // edge tangents (reshape scratch)
541  // Kokkos::DynRankView<scratchValueType,scratchProperties...> edgeTangents(scratch.data(),
542  // inputJac.extent(0),
543  // inputJac.extent(1),
544  // inputJac.extent(2));
545  auto vcprop = Kokkos::common_view_alloc_prop(scratch);
546  //typedef Kokkos::DynRankView<scratchValueType, typename decltype(scratch)::memory_space> viewType;
547  typedef Kokkos::DynRankView<scratchValueType, SpT> viewType;
548  viewType edgeTangents(Kokkos::view_wrap(scratch.data(), vcprop),
549  inputJac.extent(0),
550  inputJac.extent(1),
551  inputJac.extent(2));
552 
553  // compute normals
554  CellTools<SpT>::getPhysicalEdgeTangents(edgeTangents, inputJac, whichEdge, parentCell);
555 
556  // compute lenghts of tangents
557  RealSpaceTools<SpT>::vectorNorm(outputVals, edgeTangents, NORM_TWO);
558 
559  // multiply with weights
560  ArrayTools<SpT>::scalarMultiplyDataData(outputVals, outputVals, inputWeights);
561  }
562 
563  // ------------------------------------------------------------------------------------
564 
565  template<typename SpT>
566  template<typename outputValValueType, class ...outputValProperties,
567  typename inputMeasureValueType, class ...inputMeasureProperties,
568  typename inputValValueType, class ...inputValProperteis>
569  void
571  multiplyMeasure( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
572  const Kokkos::DynRankView<inputMeasureValueType,inputMeasureProperties...> inputMeasure,
573  const Kokkos::DynRankView<inputValValueType, inputValProperteis...> inputVals ) {
574  scalarMultiplyDataField( outputVals,
575  inputMeasure,
576  inputVals );
577  }
578 
579  // ------------------------------------------------------------------------------------
580 
581  template<typename SpT>
582  template<typename outputFieldValueType, class ...outputFieldProperties,
583  typename inputDataValueType, class ...inputDataProperties,
584  typename inputFieldValueType, class ...inputFieldProperties>
585  void
587  scalarMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
588  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
589  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields,
590  const bool reciprocal ) {
592  inputData,
593  inputFields,
594  reciprocal );
595  }
596 
597  // ------------------------------------------------------------------------------------
598 
599  template<typename SpT>
600  template<typename outputDataValuetype, class ...outputDataProperties,
601  typename inputDataLeftValueType, class ...inputDataLeftProperties,
602  typename inputDataRightValueType, class ...inputDataRightProperties>
603  void
605  scalarMultiplyDataData( Kokkos::DynRankView<outputDataValuetype, outputDataProperties...> outputData,
606  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
607  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight,
608  const bool reciprocal ) {
610  inputDataLeft,
611  inputDataRight,
612  reciprocal );
613  }
614 
615  // ------------------------------------------------------------------------------------
616 
617  template<typename SpT>
618  template<typename outputFieldValueType, class ...outputFieldProperties,
619  typename inputDataValueType, class ...inputDataProperties,
620  typename inputFieldValueType, class ...inputFieldProperties>
621  void
623  dotMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
624  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
625  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields ) {
627  inputData,
628  inputFields );
629  }
630 
631  // ------------------------------------------------------------------------------------
632 
633  template<typename SpT>
634  template<typename outputDataValueType, class ...outputDataProperties,
635  typename inputDataLeftValueType, class ...inputDataLeftProperties,
636  typename inputDataRightValueType, class ...inputDataRightProperties>
637  void
639  dotMultiplyDataData( Kokkos::DynRankView<outputDataValueType, outputDataProperties...> outputData,
640  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
641  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight ) {
643  inputDataLeft,
644  inputDataRight );
645  }
646 
647  // ------------------------------------------------------------------------------------
648 
649  template<typename SpT>
650  template<typename outputFieldValueType, class ...outputFieldProperties,
651  typename inputDataValueType, class ...inputDataProperties,
652  typename inputFieldValueType, class ...inputFieldProperties>
653  void
655  vectorMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
656  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
657  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields ) {
658  const auto outRank = outputFields.rank();
659  switch (outRank) {
660  case 3:
661  case 4:
663  inputData,
664  inputFields );
665  break;
666  case 5:
668  inputData,
669  inputFields );
670  break;
671  default: {
672  INTREPID2_TEST_FOR_EXCEPTION( outRank < 3 && outRank > 5, std::runtime_error,
673  ">>> ERROR (FunctionSpaceTools::vectorMultiplyDataField): Output container must have rank 3, 4 or 5.");
674  }
675  }
676  }
677 
678  // ------------------------------------------------------------------------------------
679 
680  template<typename SpT>
681  template<typename outputDataValueType, class ...outputDataProperties,
682  typename inputDataLeftValueType, class ...inputDataLeftProperties,
683  typename inputDataRightValueType, class ...inputDataRightProperties>
684  void
686  vectorMultiplyDataData( Kokkos::DynRankView<outputDataValueType, outputDataProperties...> outputData,
687  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
688  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight ) {
689  const auto outRank = outputData.rank();
690  switch (outRank) {
691  case 2:
692  case 3:
694  inputDataLeft,
695  inputDataRight );
696  break;
697  case 4:
699  inputDataLeft,
700  inputDataRight );
701  break;
702  default: {
703  INTREPID2_TEST_FOR_EXCEPTION( outRank < 2 && outRank > 4, std::runtime_error,
704  ">>> ERROR (FunctionSpaceTools::vectorMultiplyDataData): Output container must have rank 2, 3 or 4.");
705  }
706  }
707  }
708 
709  // ------------------------------------------------------------------------------------
710 
711  template<typename SpT>
712  template<typename outputFieldValueType, class ...outputFieldProperties,
713  typename inputDataValueType, class ...inputDataProperties,
714  typename inputFieldValueType, class ...inputFieldProperties>
715  void
717  tensorMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
718  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
719  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields,
720  const char transpose ) {
721 
722  const auto outRank = outputFields.rank();
723  switch (outRank) {
724  case 4:
726  inputData,
727  inputFields,
728  transpose );
729  break;
730  case 5:
732  inputData,
733  inputFields,
734  transpose );
735  break;
736  default: {
737  INTREPID2_TEST_FOR_EXCEPTION( outRank < 4 && outRank > 5, std::runtime_error,
738  ">>> ERROR (FunctionSpaceTools::tensorMultiplyDataField): Output container must have rank 4 or 5.");
739  }
740  }
741  }
742 
743  // ------------------------------------------------------------------------------------
744 
745  template<typename SpT>
746  template<typename outputDataValueType, class ...outputDataProperties,
747  typename inputDataLeftValueType, class ...inputDataLeftProperties,
748  typename inputDataRightValueType, class ...inputDataRightProperties>
749  void
751  tensorMultiplyDataData( Kokkos::DynRankView<outputDataValueType, outputDataProperties...> outputData,
752  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
753  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight,
754  const char transpose ) {
755  const auto outRank = outputData.rank();
756  switch (outRank) {
757  case 3:
759  inputDataLeft,
760  inputDataRight,
761  transpose );
762  break;
763  case 4:
765  inputDataLeft,
766  inputDataRight,
767  transpose );
768  break;
769  default: {
770  INTREPID2_TEST_FOR_EXCEPTION( outRank < 4 && outRank > 5, std::runtime_error,
771  ">>> ERROR (FunctionSpaceTools::tensorMultiplyDataField): Output container must have rank 4 or 5.");
772  }
773  }
774  }
775 
776  // ------------------------------------------------------------------------------------
777 
778  namespace FunctorFunctionSpaceTools {
779 
783  template<typename inoutOperatorViewType,
784  typename fieldSignViewType>
786  inoutOperatorViewType _inoutOperator;
787  const fieldSignViewType _fieldSigns;
788 
789  KOKKOS_INLINE_FUNCTION
790  F_applyLeftFieldSigns( inoutOperatorViewType inoutOperator_,
791  fieldSignViewType fieldSigns_ )
792  : _inoutOperator(inoutOperator_), _fieldSigns(fieldSigns_) {}
793 
794  KOKKOS_INLINE_FUNCTION
795  void operator()(const ordinal_type cl) const {
796  const ordinal_type nlbf = _inoutOperator.extent(1);
797  const ordinal_type nrbf = _inoutOperator.extent(2);
798 
799  for (ordinal_type lbf=0;lbf<nlbf;++lbf)
800  for (ordinal_type rbf=0;rbf<nrbf;++rbf)
801  _inoutOperator(cl, lbf, rbf) *= _fieldSigns(cl, lbf);
802  }
803  };
804  }
805 
806  template<typename SpT>
807  template<typename inoutOperatorValueType, class ...inoutOperatorProperties,
808  typename fieldSignValueType, class ...fieldSignProperties>
809  void
811  applyLeftFieldSigns( Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...> inoutOperator,
812  const Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSigns ) {
813 
814 #ifdef HAVE_INTREPID2_DEBUG
815  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.rank() != 3, std::invalid_argument,
816  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Input operator container must have rank 3.");
817  INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
818  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Input field signs container must have rank 2.");
819  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(0) != fieldSigns.extent(0), std::invalid_argument,
820  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Zeroth dimensions (number of cells) of the operator and field signs containers must agree!");
821  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(1) != fieldSigns.extent(1), std::invalid_argument,
822  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): First dimensions (number of left fields) of the operator and field signs containers must agree!");
823 #endif
824 
825  typedef Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...> inoutOperatorViewType;
826  typedef Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSignViewType;
828  <inoutOperatorViewType,fieldSignViewType> FunctorType;
829  typedef typename ExecSpace<typename inoutOperatorViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
830 
831  const ordinal_type C = inoutOperator.extent(0);
832  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
833  Kokkos::parallel_for( policy, FunctorType(inoutOperator, fieldSigns) );
834  }
835 
836  // ------------------------------------------------------------------------------------
837 
838  namespace FunctorFunctionSpaceTools {
842  template<typename inoutOperatorViewType,
843  typename fieldSignViewType>
845  inoutOperatorViewType _inoutOperator;
846  const fieldSignViewType _fieldSigns;
847 
848  KOKKOS_INLINE_FUNCTION
849  F_applyRightFieldSigns( inoutOperatorViewType inoutOperator_,
850  fieldSignViewType fieldSigns_ )
851  : _inoutOperator(inoutOperator_), _fieldSigns(fieldSigns_) {}
852 
853  KOKKOS_INLINE_FUNCTION
854  void operator()(const ordinal_type cl) const {
855  const ordinal_type nlbf = _inoutOperator.extent(1);
856  const ordinal_type nrbf = _inoutOperator.extent(2);
857 
858  for (ordinal_type lbf=0;lbf<nlbf;++lbf)
859  for (ordinal_type rbf=0;rbf<nrbf;++rbf)
860  _inoutOperator(cl, lbf, rbf) *= _fieldSigns(cl, rbf);
861  }
862  };
863  }
864 
865  template<typename SpT>
866  template<typename inoutOperatorValueType, class ...inoutOperatorProperties,
867  typename fieldSignValueType, class ...fieldSignProperties>
868  void
870  applyRightFieldSigns( Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...> inoutOperator,
871  const Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSigns ) {
872 
873 #ifdef HAVE_INTREPID2_DEBUG
874  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.rank() != 3, std::invalid_argument,
875  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Input operator container must have rank 3.");
876  INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
877  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Input field signs container must have rank 2.");
878  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(0) != fieldSigns.extent(0), std::invalid_argument,
879  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Zeroth dimensions (number of cells) of the operator and field signs containers must agree!");
880  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(2) != fieldSigns.extent(1), std::invalid_argument,
881  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Second dimension of the operator container and first dimension of the field signs container (number of right fields) must agree!");
882 #endif
883 
884  typedef Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...> inoutOperatorViewType;
885  typedef Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSignViewType;
887  <inoutOperatorViewType,fieldSignViewType> FunctorType;
888  typedef typename ExecSpace<typename inoutOperatorViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
889 
890  const ordinal_type C = inoutOperator.extent(0);
891  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
892  Kokkos::parallel_for( policy, FunctorType(inoutOperator, fieldSigns) );
893  }
894 
895  // ------------------------------------------------------------------------------------
896 
897  namespace FunctorFunctionSpaceTools {
901  template<typename inoutFunctionViewType,
902  typename fieldSignViewType>
904  inoutFunctionViewType _inoutFunction;
905  const fieldSignViewType _fieldSigns;
906 
907  KOKKOS_INLINE_FUNCTION
908  F_applyFieldSigns( inoutFunctionViewType inoutFunction_,
909  fieldSignViewType fieldSigns_)
910  : _inoutFunction(inoutFunction_), _fieldSigns(fieldSigns_) {}
911 
912  KOKKOS_INLINE_FUNCTION
913  void operator()(const ordinal_type cl) const {
914  const ordinal_type nbfs = _inoutFunction.extent(1);
915  const ordinal_type npts = _inoutFunction.extent(2);
916  const ordinal_type iend = _inoutFunction.extent(3);
917  const ordinal_type jend = _inoutFunction.extent(4);
918 
919  for (ordinal_type bf=0;bf<nbfs;++bf)
920  for (ordinal_type pt=0;pt<npts;++pt)
921  for (ordinal_type i=0;i<iend;++i)
922  for (ordinal_type j=0;j<jend;++j)
923  _inoutFunction(cl, bf, pt, i, j) *= _fieldSigns(cl, bf);
924  }
925  };
926  }
927 
928  template<typename SpT>
929  template<typename inoutFunctionValueType, class ...inoutFunctionProperties,
930  typename fieldSignValueType, class ...fieldSignProperties>
931  void
933  applyFieldSigns( Kokkos::DynRankView<inoutFunctionValueType,inoutFunctionProperties...> inoutFunction,
934  const Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSigns ) {
935 
936 #ifdef HAVE_INTREPID2_DEBUG
937  INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.rank() < 2 || inoutFunction.rank() > 5, std::invalid_argument,
938  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Input function container must have rank 2, 3, 4, or 5.");
939  INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
940  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Input field signs container must have rank 2.");
941  INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.extent(0) != fieldSigns.extent(0), std::invalid_argument,
942  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Zeroth dimensions (number of integration domains) of the function and field signs containers must agree!");
943  INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.extent(1) != fieldSigns.extent(1), std::invalid_argument,
944  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): First dimensions (number of fields) of the function and field signs containers must agree!");
945 
946 #endif
947 
948  typedef Kokkos::DynRankView<inoutFunctionValueType,inoutFunctionProperties...> inoutFunctionViewType;
949  typedef Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSignViewType;
951  <inoutFunctionViewType,fieldSignViewType> FunctorType;
952  typedef typename ExecSpace<typename inoutFunctionViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
953 
954  const ordinal_type C = inoutFunction.extent(0);
955  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
956  Kokkos::parallel_for( policy, FunctorType(inoutFunction, fieldSigns) );
957  }
958 
959  // ------------------------------------------------------------------------------------
960 
961  namespace FunctorFunctionSpaceTools {
962 
966  template<typename outputPointViewType,
967  typename inputCoeffViewType,
968  typename inputFieldViewType>
969  struct F_evaluate {
970  outputPointViewType _outputPointVals;
971  const inputCoeffViewType _inputCoeffs;
972  const inputFieldViewType _inputFields;
973 
974  KOKKOS_INLINE_FUNCTION
975  F_evaluate( outputPointViewType outputPointVals_,
976  inputCoeffViewType inputCoeffs_,
977  inputFieldViewType inputFields_ )
978  : _outputPointVals(outputPointVals_), _inputCoeffs(inputCoeffs_), _inputFields(inputFields_) {}
979 
980  KOKKOS_INLINE_FUNCTION
981  void operator()(const ordinal_type cl) const {
982  const ordinal_type nbfs = _inputFields.extent(1);
983  const ordinal_type npts = _inputFields.extent(2);
984 
985  const ordinal_type iend = _inputFields.extent(3);
986  const ordinal_type jend = _inputFields.extent(4);
987 
988  for (ordinal_type bf=0;bf<nbfs;++bf)
989  for (ordinal_type pt=0;pt<npts;++pt)
990  for (ordinal_type i=0;i<iend;++i)
991  for (ordinal_type j=0;j<jend;++j)
992  _outputPointVals(cl, pt, i, j) += _inputCoeffs(cl, bf) * _inputFields(cl, bf, pt, i, j);
993  }
994  };
995  }
996 
997  template<typename SpT>
998  template<typename outputPointValueType, class ...outputPointProperties,
999  typename inputCoeffValueType, class ...inputCoeffProperties,
1000  typename inputFieldValueType, class ...inputFieldProperties>
1001  void
1003  evaluate( Kokkos::DynRankView<outputPointValueType,outputPointProperties...> outputPointVals,
1004  const Kokkos::DynRankView<inputCoeffValueType, inputCoeffProperties...> inputCoeffs,
1005  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields ) {
1006 
1007 #ifdef HAVE_INTREPID2_DEBUG
1008  INTREPID2_TEST_FOR_EXCEPTION( inputFields.rank() < 3 || inputFields.rank() > 5, std::invalid_argument,
1009  ">>> ERROR (FunctionSpaceTools::evaluate): Input fields container must have rank 3, 4, or 5.");
1010  INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.rank() != 2, std::invalid_argument,
1011  ">>> ERROR (FunctionSpaceTools::evaluate): Input coefficient container must have rank 2.");
1012  INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.rank() != (inputFields.rank()-1), std::invalid_argument,
1013  ">>> ERROR (FunctionSpaceTools::evaluate): Output values container must have rank one less than the rank of the input fields container.");
1014  INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.extent(0) != inputFields.extent(0), std::invalid_argument,
1015  ">>> ERROR (FunctionSpaceTools::evaluate): Zeroth dimensions (number of cells) of the coefficient and fields input containers must agree!");
1016  INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.extent(1) != inputFields.extent(1), std::invalid_argument,
1017  ">>> ERROR (FunctionSpaceTools::evaluate): First dimensions (number of fields) of the coefficient and fields input containers must agree!");
1018  INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.extent(0) != inputFields.extent(0), std::invalid_argument,
1019  ">>> ERROR (FunctionSpaceTools::evaluate): Zeroth dimensions (number of cells) of the input fields container and the output values container must agree!");
1020  for (size_type i=1;i<outputPointVals.rank();++i)
1021  INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.extent(i) != inputFields.extent(i+1), std::invalid_argument,
1022  ">>> ERROR (FunctionSpaceTools::evaluate): outputPointVals dimension(i) does not match to inputFields dimension(i+1).");
1023 #endif
1024 
1025  typedef Kokkos::DynRankView<outputPointValueType,outputPointProperties...> outputPointValViewType;
1026  typedef Kokkos::DynRankView<inputCoeffValueType, inputCoeffProperties...> inputCoeffViewType;
1027  typedef Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFieldViewType;
1029  <outputPointValViewType,inputCoeffViewType,inputFieldViewType> FunctorType;
1030  typedef typename ExecSpace<typename inputCoeffViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
1031 
1032  const ordinal_type C = inputFields.extent(0);
1033  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
1034  Kokkos::parallel_for( policy, FunctorType(outputPointVals, inputCoeffs, inputFields) );
1035  }
1036 
1037  // ------------------------------------------------------------------------------------
1038 
1039 } // end namespace Intrepid2
1040 
1041 #endif
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 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 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...
Functor for calculation HGRADtransformGRAD, see Intrepid2::FunctionSpaceTools for more...
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...
Functor for applyLeftFieldSigns, see Intrepid2::FunctionSpaceTools for more.
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...
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 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 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 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 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 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 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 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 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...
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 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 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 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...
Functor for applyRightFieldSigns, see Intrepid2::FunctionSpaceTools for more.
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...
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 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 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 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 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 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...
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 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 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 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 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 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 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 HGRADtransformGRAD(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< jacobianInverseValueType, jacobianInverseProperties...> jacobianInverse, const Kokkos::DynRankView< inputValValueType, inputValProperties...> inputVals)
Transformation of a gradient field in the H-grad space, defined at points on a reference cell...
static void clone(Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input)
Clone input array.
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...
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 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 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 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 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 ...
Functor to evaluate functions, see Intrepid2::FunctionSpaceTools for more.
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 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 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 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...
Functor for calculation of cell measure, see Intrepid2::FunctionSpaceTools for more.
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 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 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 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.
Functor for applyFieldSigns, see Intrepid2::FunctionSpaceTools for more.
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 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 bool computeCellMeasure(Kokkos::DynRankView< outputValValueType, outputValProperties...> outputVals, const Kokkos::DynRankView< inputDetValueType, inputDetPropertes...> inputDet, const Kokkos::DynRankView< inputWeightValueType, inputWeightPropertes...> inputWeights)
Returns the weighted integration measures outputVals with dimensions (C,P) used for the computation o...