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 
54 
55 #include "Teuchos_TimeMonitor.hpp"
56 
57 namespace Intrepid2 {
58  // ------------------------------------------------------------------------------------
59  template<typename DeviceType>
60  template<typename outputValueType, class ...outputProperties,
61  typename inputValueType, class ...inputProperties>
62  void
64  HGRADtransformVALUE( Kokkos::DynRankView<outputValueType,outputProperties...> output,
65  const Kokkos::DynRankView<inputValueType, inputProperties...> input ) {
66  if(output.rank() == input.rank()) {
67 #ifdef HAVE_INTREPID2_DEBUG
68  {
69  for (size_type i=0;i< input.rank();++i) {
70  INTREPID2_TEST_FOR_EXCEPTION( (input.extent(i) != output.extent(i)), std::invalid_argument,
71  ">>> ERROR (FunctionSpaceTools::HGRADtransformVALUE): Dimensions of input and output fields containers do not match.");
72  }
73  }
74 #endif
75  RealSpaceTools<DeviceType>::clone(output, input);
76  }
77  else
79  }
80 
81  template<typename DeviceType>
82  template<typename outputValueType, class ...outputProperties,
83  typename inputValueType, class ...inputProperties>
84  void
86  mapHGradDataFromPhysToRef( Kokkos::DynRankView<outputValueType,outputProperties...> output,
87  const Kokkos::DynRankView<inputValueType, inputProperties...> input ) {
88  if(output.rank() == input.rank()) {
89 #ifdef HAVE_INTREPID2_DEBUG
90  {
91  for (size_type i=0;i< input.rank();++i) {
92  INTREPID2_TEST_FOR_EXCEPTION( (input.extent(i) != output.extent(i)), std::invalid_argument,
93  ">>> ERROR (FunctionSpaceTools::mapHGradDataFromPhysToRef): Dimensions of input and output fields containers do not match.");
94  }
95  }
96 #endif
97  RealSpaceTools<DeviceType>::clone(output, input);
98  }
99  else
101  }
102 
103  template<typename DeviceType>
104  template<typename outputValueType, class ...outputProperties,
105  typename inputValueType, class ...inputProperties>
106  void
108  mapHGradDataFromPhysSideToRefSide( Kokkos::DynRankView<outputValueType,outputProperties...> output,
109  const Kokkos::DynRankView<inputValueType, inputProperties...> input ) {
110  if(output.rank() == input.rank()) {
111 #ifdef HAVE_INTREPID2_DEBUG
112  {
113  for (size_type i=0;i< input.rank();++i) {
114  INTREPID2_TEST_FOR_EXCEPTION( (input.extent(i) != output.extent(i)), std::invalid_argument,
115  ">>> ERROR (FunctionSpaceTools::mapHGradDataSideFromPhysToRefSide): Dimensions of input and output fields containers do not match.");
116  }
117  }
118 #endif
119  RealSpaceTools<DeviceType>::clone(output, input);
120  }
121  else
123  }
124 
125  // ------------------------------------------------------------------------------------
126 
127  namespace FunctorFunctionSpaceTools {
131  }
132 
133  template<typename DeviceType>
134  template<typename OutputValViewType,
135  typename JacobianInverseViewType,
136  typename InputValViewType>
137  void
139  HGRADtransformGRAD( OutputValViewType outputVals,
140  const JacobianInverseViewType jacobianInverse,
141  const InputValViewType inputVals ) {
142  return HCURLtransformVALUE(outputVals, jacobianInverse, inputVals);
143  }
144 
145  // ------------------------------------------------------------------------------------
146 
147  template<typename DeviceType>
148  template<typename outputValValueType, class ...outputValProperties,
149  typename jacobianInverseValueType, class ...jacobianInverseProperties,
150  typename inputValValueType, class ...inputValProperties>
151  void
153  HCURLtransformVALUE( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
154  const Kokkos::DynRankView<jacobianInverseValueType,jacobianInverseProperties...> jacobianInverse,
155  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
156  ArrayTools<DeviceType>::matvecProductDataField(outputVals, jacobianInverse, inputVals, 'T');
157  }
158 
159  template<typename DeviceType>
160  template<typename outputValValueType, class ...outputValProperties,
161  typename jacobianValueType, class ...jacobianProperties,
162  typename inputValValueType, class ...inputValProperties>
163  void
165  mapHCurlDataFromPhysToRef( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
166  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
167  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
168  ArrayTools<DeviceType>::matvecProductDataData(outputVals, jacobian, inputVals, 'T');
169  }
170 
171 
172  namespace FunctorFunctionSpaceTools {
173 
174  template<typename outViewType,
175  typename inputViewType,
176  typename metricViewType
177  >
179  outViewType output_;
180  const inputViewType input_;
181  const metricViewType metricTensorDet_;
182 
183  KOKKOS_INLINE_FUNCTION
184  F_negativeWeighted2dInputCrossK( outViewType output,
185  const inputViewType input,
186  const metricViewType metricTensorDet)
187  : output_(output), input_(input), metricTensorDet_(metricTensorDet){};
188 
189  KOKKOS_INLINE_FUNCTION
190  void operator()(const size_type ic) const {
191  for (size_t pt=0; pt < input_.extent(1); pt++) {
192  auto measure = std::sqrt(metricTensorDet_(ic,pt));
193  output_(ic,pt,0) = - measure*input_(ic,pt,1);
194  output_(ic,pt,1) = measure*input_(ic,pt,0);
195  }
196  }
197  };
198  }
199 
200  template<typename DeviceType>
201  template<typename outputValValueType, class ...outputValProperties,
202  typename tangentsValueType, class ...tangentsProperties,
203  typename metricTensorInvValueType, class ...metricTensorInvProperties,
204  typename metricTensorDetValueType, class ...metricTensorDetProperties,
205  typename inputValValueType, class ...inputValProperties>
206  void
208  mapHCurlDataCrossNormalFromPhysSideToRefSide( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
209  const Kokkos::DynRankView<tangentsValueType, tangentsProperties...> tangents,
210  const Kokkos::DynRankView<metricTensorInvValueType,metricTensorInvProperties...> metricTensorInv,
211  const Kokkos::DynRankView<metricTensorDetValueType,metricTensorDetProperties...> metricTensorDet,
212  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
213  auto work = Kokkos::DynRankView<outputValValueType, outputValProperties...>("work", outputVals.extent(0), outputVals.extent(1),outputVals.extent(2));
214  ArrayTools<DeviceType>::matvecProductDataData(outputVals, tangents, inputVals, 'T');
215  typename DeviceType::execution_space().fence();
216  ArrayTools<DeviceType>::matvecProductDataData(work, metricTensorInv, outputVals);
217  typename DeviceType::execution_space().fence();
219  Kokkos::parallel_for(Kokkos::RangePolicy<ExecSpaceType>(0,outputVals.extent(0)), FunctorType(outputVals, work, metricTensorDet) );
220  }
221 
222 
223  namespace FunctorFunctionSpaceTools {
224 
225  template<typename outViewType,
226  typename inputViewType,
227  typename metricViewType
228  >
229  struct F_weighedInput {
230  outViewType output_;
231  const inputViewType input_;
232  const metricViewType metricTensorDet_;
233  const double scaling_;
234 
235  KOKKOS_INLINE_FUNCTION
236  F_weighedInput( outViewType output,
237  const inputViewType input,
238  const metricViewType metricTensorDet,
239  const double scaling)
240  : output_(output), input_(input), metricTensorDet_(metricTensorDet), scaling_(scaling){};
241 
242  KOKKOS_INLINE_FUNCTION
243  void operator()(const size_type ic) const {
244  for (size_t pt=0; pt < input_.extent(1); pt++) {
245  auto measure = std::sqrt(metricTensorDet_(ic,pt));
246  output_.access(ic,pt) = scaling_ * measure * input_.access(ic,pt);
247  }
248  }
249  };
250  }
251 
252  template<typename DeviceType>
253  template<typename outputValValueType, class ...outputValProperties,
254  typename jacobianDetValueType, class ...jacobianDetProperties,
255  typename inputValValueType, class ...inputValProperties>
256  void
259  Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
260  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> metricTensorDet,
261  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
263  Kokkos::parallel_for(Kokkos::RangePolicy<ExecSpaceType>(0,outputVals.extent(0)), FunctorType(outputVals, inputVals, metricTensorDet, -1.0) );
264  }
265 
266  // ------------------------------------------------------------------------------------
267 
268  template<typename DeviceType>
269  template<typename outputValValueType, class ...outputValProperties,
270  typename jacobianValueType, class ...jacobianProperties,
271  typename jacobianDetValueType, class ...jacobianDetProperties,
272  typename inputValValueType, class ...inputValProperties>
273  void
275  HCURLtransformCURL( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
276  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
277  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
278  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
279  if(jacobian.data()==NULL || jacobian.extent(2)==2) //2D case
280  return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
281  else
282  return HDIVtransformVALUE(outputVals, jacobian, jacobianDet, inputVals);
283  }
284 
285  // ------------------------------------------------------------------------------------
286 
287  template<typename DeviceType>
288  template<typename outputValValueType, class ...outputValProperties,
289  typename jacobianDetValueType, class ...jacobianDetProperties,
290  typename inputValValueType, class ...inputValProperties>
291  void
293  HCURLtransformCURL( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
294  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
295  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
296 #ifdef HAVE_INTREPID2_DEBUG
297  {
298  INTREPID2_TEST_FOR_EXCEPTION( outputVals.rank() == 4, std::invalid_argument,
299  ">>> ERROR (FunctionSpaceTools::HCURLtransformCURL): Output rank must have rank 3.\n If these are 3D fields, then use the appropriate overload of this function.");
300  }
301 #endif
302  return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
303  }
304 
305  // ------------------------------------------------------------------------------------
306 
307  template<typename DeviceType>
308  template<typename outputValValueType, class ...outputValProperties,
309  typename jacobianValueType, class ...jacobianProperties,
310  typename jacobianDetValueType, class ...jacobianDetProperties,
311  typename inputValValueType, class ...inputValProperties>
312  void
314  HGRADtransformCURL( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
315  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
316  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
317  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
318 #ifdef HAVE_INTREPID2_DEBUG
319  {
320  INTREPID2_TEST_FOR_EXCEPTION( outputVals.extent(3)!=2, std::invalid_argument,
321  ">>> ERROR (FunctionSpaceTools::HGRADtransformCURL):\n output field is 3D by the function is meant for 2D fields");
322  }
323 #endif
324  return HDIVtransformVALUE(outputVals, jacobian, jacobianDet, inputVals);
325  }
326 
327  // ------------------------------------------------------------------------------------
328 
329  template<typename DeviceType>
330  template<typename outputValValueType, class ...outputValProperties,
331  typename jacobianValueType, class ...jacobianProperties,
332  typename jacobianDetValueType, class ...jacobianDetProperties,
333  typename inputValValueType, class ...inputValProperties>
334  void
336  HDIVtransformVALUE( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
337  const Kokkos::DynRankView<jacobianValueType, jacobianProperties...> jacobian,
338  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
339  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
340  ArrayTools<DeviceType>::matvecProductDataField(outputVals, jacobian, inputVals, 'N');
341  ArrayTools<DeviceType>::scalarMultiplyDataField(outputVals, jacobianDet, outputVals, true);
342  }
343 
344  template<typename DeviceType>
345  template<typename outputValValueType, class ...outputValProperties,
346  typename jacobianInverseValueType, class ...jacobianInverseProperties,
347  typename jacobianDetValueType, class ...jacobianDetProperties,
348  typename inputValValueType, class ...inputValProperties>
349  void
351  mapHDivDataFromPhysToRef( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
352  const Kokkos::DynRankView<jacobianInverseValueType, jacobianInverseProperties...> jacobianInv,
353  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
354  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
355  ArrayTools<DeviceType>::matvecProductDataData(outputVals, jacobianInv, inputVals, 'N');
356  typename DeviceType::execution_space().fence();
357  ArrayTools<DeviceType>::scalarMultiplyDataData(outputVals, jacobianDet, outputVals, false);
358  }
359 
360 
361 
362  template<typename DeviceType>
363  template<typename outputValValueType, class ...outputValProperties,
364  typename jacobianDetValueType, class ...jacobianDetProperties,
365  typename inputValValueType, class ...inputValProperties>
366  void
369  Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
370  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> metricTensorDet,
371  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
373  Kokkos::parallel_for(Kokkos::RangePolicy<ExecSpaceType>(0,outputVals.extent(0)), FunctorType(outputVals, inputVals, metricTensorDet, 1.0) );
374  }
375 
376  // ------------------------------------------------------------------------------------
377 
378  template<typename DeviceType>
379  template<typename outputValValueType, class ...outputValProperties,
380  typename jacobianDetValueType, class ...jacobianDetProperties,
381  typename inputValValueType, class ...inputValProperties>
382  void
384  HDIVtransformDIV( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
385  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
386  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
387  return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
388  }
389 
390  // ------------------------------------------------------------------------------------
391 
392  template<typename DeviceType>
393  template<typename outputValValueType, class ...outputValProperties,
394  typename jacobianDetValueType, class ...jacobianDetProperties,
395  typename inputValValueType, class ...inputValProperties>
396  void
398  HVOLtransformVALUE( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
399  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
400  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
401  ArrayTools<DeviceType>::scalarMultiplyDataField(outputVals, jacobianDet, inputVals, true);
402  }
403 
404  template<typename DeviceType>
405  template<typename outputValValueType, class ...outputValProperties,
406  typename jacobianDetValueType, class ...jacobianDetProperties,
407  typename inputValValueType, class ...inputValProperties>
408  void
410  mapHVolDataFromPhysToRef( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
411  const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
412  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
413  ArrayTools<DeviceType>::scalarMultiplyDataData(outputVals, jacobianDet, inputVals, false);
414  }
415 
416 
417 
418  // ------------------------------------------------------------------------------------
419 
420  template<typename DeviceType>
421  template<typename outputValueValueType, class ...outputValueProperties,
422  typename leftValueValueType, class ...leftValueProperties,
423  typename rightValueValueType, class ...rightValueProperties>
424  void
426  integrate( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
427  const Kokkos::DynRankView<leftValueValueType, leftValueProperties...> leftValues,
428  const Kokkos::DynRankView<rightValueValueType, rightValueProperties...> rightValues,
429  const bool sumInto ) {
430 
431 #ifdef HAVE_INTREPID2_DEBUG
432  {
433  INTREPID2_TEST_FOR_EXCEPTION( leftValues.rank() < 2 ||
434  leftValues.rank() > 4, std::invalid_argument,
435  ">>> ERROR (FunctionSpaceTools::integrate): Left data must have rank 2, 3 or 4.");
436  INTREPID2_TEST_FOR_EXCEPTION( outputValues.rank() < 1 ||
437  outputValues.rank() > 3, std::invalid_argument,
438  ">>> ERROR (FunctionSpaceTools::integrate): Output values must have rank 1, 2 or 3.");
439  }
440 #endif
441 
442  const ordinal_type outRank = outputValues.rank();
443  const ordinal_type leftRank = leftValues.rank();
444  const ordinal_type mode = outRank*10 + leftRank;
445 
446  switch (mode) {
447  // DataData
448  case 12:
450  leftValues,
451  rightValues,
452  sumInto );
453  break;
454  case 13:
456  leftValues,
457  rightValues,
458  sumInto );
459  break;
460  case 14:
462  leftValues,
463  rightValues,
464  sumInto );
465  break;
466 
467  // DataField
468  case 22:
470  leftValues,
471  rightValues,
472  sumInto );
473  break;
474  case 23:
476  leftValues,
477  rightValues,
478  sumInto );
479  break;
480  case 24:
482  leftValues,
483  rightValues,
484  sumInto );
485  break;
486 
487  // FieldField
488  case 33:
490  leftValues,
491  rightValues,
492  sumInto );
493  break;
494  case 34:
496  leftValues,
497  rightValues,
498  sumInto );
499  break;
500  case 35:
502  leftValues,
503  rightValues,
504  sumInto );
505  break;
506  default: {
507  INTREPID2_TEST_FOR_EXCEPTION( outRank < 1 || outRank > 3, std::runtime_error,
508  ">>> ERROR (FunctionSpaceTools::integrate): outRank must be 1,2, or 3.");
509  INTREPID2_TEST_FOR_EXCEPTION( leftRank < 2 || leftRank > 4, std::runtime_error,
510  ">>> ERROR (FunctionSpaceTools::integrate): leftRank must be 1,2, 3 or 4.");
511  }
512  }
513  }
514 
515  // ------------------------------------------------------------------------------------
516  namespace FunctorFunctionSpaceTools {
520  template<typename outputValViewType,
521  typename inputDetViewType,
522  typename inputWeightViewType>
524  outputValViewType _outputVals;
525  const inputDetViewType _inputDet;
526  const inputWeightViewType _inputWeight;
527 
528  KOKKOS_INLINE_FUNCTION
529  F_computeCellMeasure( outputValViewType outputVals_,
530  inputDetViewType inputDet_,
531  inputWeightViewType inputWeight_)
532  : _outputVals(outputVals_),
533  _inputDet(inputDet_),
534  _inputWeight(inputWeight_) {}
535 
536  typedef ordinal_type value_type;
537 
538 // KOKKOS_INLINE_FUNCTION
539 // void init( value_type &dst ) const {
540 // dst = false;
541 // }
542 
543 // KOKKOS_INLINE_FUNCTION
544 // void join( volatile value_type &dst,
545 // const volatile value_type &src ) const {
546 // dst |= src;
547 // }
548 
549  KOKKOS_INLINE_FUNCTION
550  void operator()(const size_type cl, value_type &dst) const {
551  // negative jacobian check
552  const bool hasNegativeDet = (_inputDet(cl, 0) < 0.0);
553  dst |= hasNegativeDet;
554 
555  // make it absolute
556  const auto sign = (hasNegativeDet ? -1.0 : 1.0);
557  const ordinal_type pt_end = _outputVals.extent(1);
558  for (ordinal_type pt=0;pt<pt_end;++pt)
559  _outputVals(cl, pt) = sign*_inputDet(cl, pt)*_inputWeight(pt);
560  }
561  };
562  }
563 
564  template<typename DeviceType>
565  template<typename OutputValViewType,
566  typename InputDetViewType,
567  typename InputWeightViewType>
568  bool
570  computeCellMeasure( OutputValViewType outputVals,
571  const InputDetViewType inputDet,
572  const InputWeightViewType inputWeights ) {
573 #ifdef HAVE_INTREPID2_DEBUG
574  {
575  INTREPID2_TEST_FOR_EXCEPTION( rank(inputDet) != 2 ||
576  rank(inputWeights) != 1 ||
577  rank(outputVals) != 2, std::invalid_argument,
578  ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Ranks are not compatible.");
579  INTREPID2_TEST_FOR_EXCEPTION( outputVals.extent(0) != inputDet.extent(0), std::invalid_argument,
580  ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Cell dimension does not match.");
581  INTREPID2_TEST_FOR_EXCEPTION( inputDet.extent(1) != outputVals.extent(1) ||
582  inputWeights.extent(0) != outputVals.extent(1), std::invalid_argument,
583  ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Point dimension does not match.");
584  }
585 #endif
586  constexpr bool are_accessible =
587  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
588  typename decltype(outputVals)::memory_space>::accessible &&
589  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
590  typename decltype(inputDet)::memory_space>::accessible &&
591  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
592  typename decltype(inputWeights)::memory_space>::accessible;
593  static_assert(are_accessible, "FunctionSpaceTools<DeviceType>::computeCellMeasure(..): input/output views' memory spaces are not compatible with DeviceType");
594 
596  <decltype(outputVals),decltype(inputDet),decltype(inputWeights)>;
597 
598  const ordinal_type C = inputDet.extent(0);
599  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
600 
601  typename FunctorType::value_type hasNegativeDet = false;
602  Kokkos::parallel_reduce( policy, FunctorType(outputVals, inputDet, inputWeights), hasNegativeDet );
603 
604  return hasNegativeDet;
605  }
606 
607  // ------------------------------------------------------------------------------------
608 
609  template<typename DeviceType>
610  template<typename outputValValueType, class ...outputValProperties,
611  typename inputJacValueType, class ...inputJacProperties,
612  typename inputWeightValueType, class ...inputWeightProperties,
613  typename scratchValueType, class ...scratchProperties>
614  void
616  computeFaceMeasure( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
617  const Kokkos::DynRankView<inputJacValueType, inputJacProperties...> inputJac,
618  const Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeights,
619  const ordinal_type whichFace,
620  const shards::CellTopology parentCell,
621  const Kokkos::DynRankView<scratchValueType, scratchProperties...> scratch ) {
622 #ifdef HAVE_INTREPID2_DEBUG
623  INTREPID2_TEST_FOR_EXCEPTION( inputJac.rank() != 4, std::invalid_argument,
624  ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Input Jacobian container must have rank 4.");
625  INTREPID2_TEST_FOR_EXCEPTION( scratch.rank() != 1, std::invalid_argument,
626  ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Scratch view imust have rank 1.");
627  INTREPID2_TEST_FOR_EXCEPTION( scratch.size() < inputJac.size(), std::invalid_argument,
628  ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Scratch storage must be greater than or equal to inputJac's one.");
629 #endif
630 
631  // face normals (reshape scratch)
632  // Kokkos::DynRankView<scratchValueType,scratchProperties...> faceNormals(scratch.data(),
633  // inputJac.extent(0),
634  // inputJac.extent(1),
635  // inputJac.extent(2));
636  auto vcprop = Kokkos::common_view_alloc_prop(scratch);
637  //typedef Kokkos::DynRankView<scratchValueType, typename decltype(scratch)::memory_space> viewType;
638  typedef Kokkos::DynRankView<scratchValueType, DeviceType> viewType;
639  viewType faceNormals(Kokkos::view_wrap(scratch.data(), vcprop),
640  inputJac.extent(0),
641  inputJac.extent(1),
642  inputJac.extent(2));
643 
644  // compute normals
645  CellTools<DeviceType>::getPhysicalFaceNormals(faceNormals, inputJac, whichFace, parentCell);
646 
647  // compute lenghts of normals
648  RealSpaceTools<DeviceType>::vectorNorm(outputVals, faceNormals, NORM_TWO);
649 
650  // multiply with weights
651  ArrayTools<DeviceType>::scalarMultiplyDataData(outputVals, outputVals, inputWeights);
652  }
653 
654  // ------------------------------------------------------------------------------------
655 
656  template<typename DeviceType>
657  template<typename outputValValueType, class ...outputValProperties,
658  typename inputJacValueType, class ...inputJacProperties,
659  typename inputWeightValueType, class ...inputWeightProperties,
660  typename scratchValueType, class ...scratchProperties>
661  void
663  computeEdgeMeasure( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
664  const Kokkos::DynRankView<inputJacValueType, inputJacProperties...> inputJac,
665  const Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeights,
666  const ordinal_type whichEdge,
667  const shards::CellTopology parentCell,
668  const Kokkos::DynRankView<scratchValueType, scratchProperties...> scratch ) {
669 #ifdef HAVE_INTREPID2_DEBUG
670  INTREPID2_TEST_FOR_EXCEPTION( (inputJac.rank() != 4), std::invalid_argument,
671  ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Input Jacobian container must have rank 4.");
672  INTREPID2_TEST_FOR_EXCEPTION( scratch.rank() != 1, std::invalid_argument,
673  ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Scratch view must have a rank 1.");
674  INTREPID2_TEST_FOR_EXCEPTION( scratch.size() < inputJac.size(), std::invalid_argument,
675  ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Scratch storage must be greater than or equal to inputJac'one.");
676 #endif
677 
678  // edge tangents (reshape scratch)
679  // Kokkos::DynRankView<scratchValueType,scratchProperties...> edgeTangents(scratch.data(),
680  // inputJac.extent(0),
681  // inputJac.extent(1),
682  // inputJac.extent(2));
683  auto vcprop = Kokkos::common_view_alloc_prop(scratch);
684  //typedef Kokkos::DynRankView<scratchValueType, typename decltype(scratch)::memory_space> viewType;
685  typedef Kokkos::DynRankView<scratchValueType, DeviceType> viewType;
686  viewType edgeTangents(Kokkos::view_wrap(scratch.data(), vcprop),
687  inputJac.extent(0),
688  inputJac.extent(1),
689  inputJac.extent(2));
690 
691  // compute normals
692  CellTools<DeviceType>::getPhysicalEdgeTangents(edgeTangents, inputJac, whichEdge, parentCell);
693 
694  // compute lenghts of tangents
695  RealSpaceTools<DeviceType>::vectorNorm(outputVals, edgeTangents, NORM_TWO);
696 
697  // multiply with weights
698  ArrayTools<DeviceType>::scalarMultiplyDataData(outputVals, outputVals, inputWeights);
699  }
700 
701  // ------------------------------------------------------------------------------------
702 
703  template<typename DeviceType>
704  template<typename outputValValueType, class ...outputValProperties,
705  typename inputMeasureValueType, class ...inputMeasureProperties,
706  typename inputValValueType, class ...inputValProperties>
707  void
709  multiplyMeasure( Kokkos::DynRankView<outputValValueType, outputValProperties...> outputVals,
710  const Kokkos::DynRankView<inputMeasureValueType,inputMeasureProperties...> inputMeasure,
711  const Kokkos::DynRankView<inputValValueType, inputValProperties...> inputVals ) {
712  scalarMultiplyDataField( outputVals,
713  inputMeasure,
714  inputVals );
715  }
716 
717  // ------------------------------------------------------------------------------------
718 
719  template<typename DeviceType>
720  template<typename outputFieldValueType, class ...outputFieldProperties,
721  typename inputDataValueType, class ...inputDataProperties,
722  typename inputFieldValueType, class ...inputFieldProperties>
723  void
725  scalarMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
726  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
727  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields,
728  const bool reciprocal ) {
730  inputData,
731  inputFields,
732  reciprocal );
733  }
734 
735  // ------------------------------------------------------------------------------------
736 
737  template<typename DeviceType>
738  template<typename outputDataValuetype, class ...outputDataProperties,
739  typename inputDataLeftValueType, class ...inputDataLeftProperties,
740  typename inputDataRightValueType, class ...inputDataRightProperties>
741  void
743  scalarMultiplyDataData( Kokkos::DynRankView<outputDataValuetype, outputDataProperties...> outputData,
744  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
745  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight,
746  const bool reciprocal ) {
748  inputDataLeft,
749  inputDataRight,
750  reciprocal );
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  dotMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
762  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
763  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields ) {
765  inputData,
766  inputFields );
767  }
768 
769  // ------------------------------------------------------------------------------------
770 
771  template<typename DeviceType>
772  template<typename outputDataValueType, class ...outputDataProperties,
773  typename inputDataLeftValueType, class ...inputDataLeftProperties,
774  typename inputDataRightValueType, class ...inputDataRightProperties>
775  void
777  dotMultiplyDataData( Kokkos::DynRankView<outputDataValueType, outputDataProperties...> outputData,
778  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
779  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight ) {
781  inputDataLeft,
782  inputDataRight );
783  }
784 
785  // ------------------------------------------------------------------------------------
786 
787  template<typename DeviceType>
788  template<typename outputFieldValueType, class ...outputFieldProperties,
789  typename inputDataValueType, class ...inputDataProperties,
790  typename inputFieldValueType, class ...inputFieldProperties>
791  void
793  vectorMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
794  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
795  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields ) {
796  const auto outRank = outputFields.rank();
797  switch (outRank) {
798  case 3:
799  case 4:
801  inputData,
802  inputFields );
803  break;
804  case 5:
806  inputData,
807  inputFields );
808  break;
809  default: {
810  INTREPID2_TEST_FOR_EXCEPTION( outRank < 3 && outRank > 5, std::runtime_error,
811  ">>> ERROR (FunctionSpaceTools::vectorMultiplyDataField): Output container must have rank 3, 4 or 5.");
812  }
813  }
814  }
815 
816  // ------------------------------------------------------------------------------------
817 
818  template<typename DeviceType>
819  template<typename outputDataValueType, class ...outputDataProperties,
820  typename inputDataLeftValueType, class ...inputDataLeftProperties,
821  typename inputDataRightValueType, class ...inputDataRightProperties>
822  void
824  vectorMultiplyDataData( Kokkos::DynRankView<outputDataValueType, outputDataProperties...> outputData,
825  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
826  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight ) {
827  const auto outRank = outputData.rank();
828  switch (outRank) {
829  case 2:
830  case 3:
832  inputDataLeft,
833  inputDataRight );
834  break;
835  case 4:
837  inputDataLeft,
838  inputDataRight );
839  break;
840  default: {
841  INTREPID2_TEST_FOR_EXCEPTION( outRank < 2 && outRank > 4, std::runtime_error,
842  ">>> ERROR (FunctionSpaceTools::vectorMultiplyDataData): Output container must have rank 2, 3 or 4.");
843  }
844  }
845  }
846 
847  // ------------------------------------------------------------------------------------
848 
849  template<typename DeviceType>
850  template<typename outputFieldValueType, class ...outputFieldProperties,
851  typename inputDataValueType, class ...inputDataProperties,
852  typename inputFieldValueType, class ...inputFieldProperties>
853  void
855  tensorMultiplyDataField( Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFields,
856  const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
857  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields,
858  const char transpose ) {
859 
860  const auto outRank = outputFields.rank();
861  switch (outRank) {
862  case 4:
864  inputData,
865  inputFields,
866  transpose );
867  break;
868  case 5:
870  inputData,
871  inputFields,
872  transpose );
873  break;
874  default: {
875  INTREPID2_TEST_FOR_EXCEPTION( outRank < 4 && outRank > 5, std::runtime_error,
876  ">>> ERROR (FunctionSpaceTools::tensorMultiplyDataField): Output container must have rank 4 or 5.");
877  }
878  }
879  }
880 
881  // ------------------------------------------------------------------------------------
882 
883  template<typename DeviceType>
884  template<typename outputDataValueType, class ...outputDataProperties,
885  typename inputDataLeftValueType, class ...inputDataLeftProperties,
886  typename inputDataRightValueType, class ...inputDataRightProperties>
887  void
889  tensorMultiplyDataData( Kokkos::DynRankView<outputDataValueType, outputDataProperties...> outputData,
890  const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
891  const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight,
892  const char transpose ) {
893  const auto outRank = outputData.rank();
894  switch (outRank) {
895  case 3:
897  inputDataLeft,
898  inputDataRight,
899  transpose );
900  break;
901  case 4:
903  inputDataLeft,
904  inputDataRight,
905  transpose );
906  break;
907  default: {
908  INTREPID2_TEST_FOR_EXCEPTION( outRank < 4 && outRank > 5, std::runtime_error,
909  ">>> ERROR (FunctionSpaceTools::tensorMultiplyDataField): Output container must have rank 4 or 5.");
910  }
911  }
912  }
913 
914  // ------------------------------------------------------------------------------------
915 
916  namespace FunctorFunctionSpaceTools {
917 
921  template<typename inoutOperatorViewType,
922  typename fieldSignViewType>
924  inoutOperatorViewType _inoutOperator;
925  const fieldSignViewType _fieldSigns;
926 
927  KOKKOS_INLINE_FUNCTION
928  F_applyLeftFieldSigns( inoutOperatorViewType inoutOperator_,
929  fieldSignViewType fieldSigns_ )
930  : _inoutOperator(inoutOperator_), _fieldSigns(fieldSigns_) {}
931 
932  KOKKOS_INLINE_FUNCTION
933  void operator()(const ordinal_type cl) const {
934  const ordinal_type nlbf = _inoutOperator.extent(1);
935  const ordinal_type nrbf = _inoutOperator.extent(2);
936 
937  for (ordinal_type lbf=0;lbf<nlbf;++lbf)
938  for (ordinal_type rbf=0;rbf<nrbf;++rbf)
939  _inoutOperator(cl, lbf, rbf) *= _fieldSigns(cl, lbf);
940  }
941  };
942  }
943 
944  template<typename DeviceType>
945  template<typename inoutOperatorValueType, class ...inoutOperatorProperties,
946  typename fieldSignValueType, class ...fieldSignProperties>
947  void
949  applyLeftFieldSigns( Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...> inoutOperator,
950  const Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSigns ) {
951 
952 #ifdef HAVE_INTREPID2_DEBUG
953  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.rank() != 3, std::invalid_argument,
954  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Input operator container must have rank 3.");
955  INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
956  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Input field signs container must have rank 2.");
957  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(0) != fieldSigns.extent(0), std::invalid_argument,
958  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Zeroth dimensions (number of cells) of the operator and field signs containers must agree!");
959  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(1) != fieldSigns.extent(1), std::invalid_argument,
960  ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): First dimensions (number of left fields) of the operator and field signs containers must agree!");
961 #endif
962 
963  constexpr bool are_accessible =
964  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
965  typename decltype(inoutOperator)::memory_space>::accessible &&
966  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
967  typename decltype(fieldSigns)::memory_space>::accessible;
968  static_assert(are_accessible, "FunctionSpaceTools<DeviceType>::applyLeftFieldSigns(..): input/output views' memory spaces are not compatible with DeviceType");
969 
971  <decltype(inoutOperator),decltype(fieldSigns)>;
972 
973  const ordinal_type C = inoutOperator.extent(0);
974  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
975  Kokkos::parallel_for( policy, FunctorType(inoutOperator, fieldSigns) );
976  }
977 
978  // ------------------------------------------------------------------------------------
979 
980  namespace FunctorFunctionSpaceTools {
984  template<typename inoutOperatorViewType,
985  typename fieldSignViewType>
987  inoutOperatorViewType _inoutOperator;
988  const fieldSignViewType _fieldSigns;
989 
990  KOKKOS_INLINE_FUNCTION
991  F_applyRightFieldSigns( inoutOperatorViewType inoutOperator_,
992  fieldSignViewType fieldSigns_ )
993  : _inoutOperator(inoutOperator_), _fieldSigns(fieldSigns_) {}
994 
995  KOKKOS_INLINE_FUNCTION
996  void operator()(const ordinal_type cl) const {
997  const ordinal_type nlbf = _inoutOperator.extent(1);
998  const ordinal_type nrbf = _inoutOperator.extent(2);
999 
1000  for (ordinal_type lbf=0;lbf<nlbf;++lbf)
1001  for (ordinal_type rbf=0;rbf<nrbf;++rbf)
1002  _inoutOperator(cl, lbf, rbf) *= _fieldSigns(cl, rbf);
1003  }
1004  };
1005  }
1006 
1007  template<typename DeviceType>
1008  template<typename inoutOperatorValueType, class ...inoutOperatorProperties,
1009  typename fieldSignValueType, class ...fieldSignProperties>
1010  void
1012  applyRightFieldSigns( Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...> inoutOperator,
1013  const Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSigns ) {
1014 
1015 #ifdef HAVE_INTREPID2_DEBUG
1016  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.rank() != 3, std::invalid_argument,
1017  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Input operator container must have rank 3.");
1018  INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
1019  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Input field signs container must have rank 2.");
1020  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(0) != fieldSigns.extent(0), std::invalid_argument,
1021  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Zeroth dimensions (number of cells) of the operator and field signs containers must agree!");
1022  INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(2) != fieldSigns.extent(1), std::invalid_argument,
1023  ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Second dimension of the operator container and first dimension of the field signs container (number of right fields) must agree!");
1024 #endif
1025 
1026  constexpr bool are_accessible =
1027  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1028  typename decltype(inoutOperator)::memory_space>::accessible &&
1029  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1030  typename decltype(fieldSigns)::memory_space>::accessible;
1031  static_assert(are_accessible, "FunctionSpaceTools<DeviceType>::applyRightFieldSigns(..): input/output views' memory spaces are not compatible with DeviceType");
1032 
1034  <decltype(inoutOperator),decltype(fieldSigns)>;
1035 
1036  const ordinal_type C = inoutOperator.extent(0);
1037  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
1038  Kokkos::parallel_for( policy, FunctorType(inoutOperator, fieldSigns) );
1039  }
1040 
1041  // ------------------------------------------------------------------------------------
1042 
1043  namespace FunctorFunctionSpaceTools {
1047  template<typename inoutFunctionViewType,
1048  typename fieldSignViewType>
1050  inoutFunctionViewType _inoutFunction;
1051  const fieldSignViewType _fieldSigns;
1052 
1053  KOKKOS_INLINE_FUNCTION
1054  F_applyFieldSigns( inoutFunctionViewType inoutFunction_,
1055  fieldSignViewType fieldSigns_)
1056  : _inoutFunction(inoutFunction_), _fieldSigns(fieldSigns_) {}
1057 
1058  KOKKOS_INLINE_FUNCTION
1059  void operator()(const ordinal_type cl) const {
1060  const ordinal_type nbfs = _inoutFunction.extent(1);
1061  const ordinal_type npts = _inoutFunction.extent(2);
1062  const ordinal_type iend = _inoutFunction.extent(3);
1063  const ordinal_type jend = _inoutFunction.extent(4);
1064 
1065  for (ordinal_type bf=0;bf<nbfs;++bf)
1066  for (ordinal_type pt=0;pt<npts;++pt)
1067  for (ordinal_type i=0;i<iend;++i)
1068  for (ordinal_type j=0;j<jend;++j)
1069  _inoutFunction(cl, bf, pt, i, j) *= _fieldSigns(cl, bf);
1070  }
1071  };
1072  }
1073 
1074  template<typename DeviceType>
1075  template<typename inoutFunctionValueType, class ...inoutFunctionProperties,
1076  typename fieldSignValueType, class ...fieldSignProperties>
1077  void
1079  applyFieldSigns( Kokkos::DynRankView<inoutFunctionValueType,inoutFunctionProperties...> inoutFunction,
1080  const Kokkos::DynRankView<fieldSignValueType, fieldSignProperties...> fieldSigns ) {
1081 
1082 #ifdef HAVE_INTREPID2_DEBUG
1083  INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.rank() < 2 || inoutFunction.rank() > 5, std::invalid_argument,
1084  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Input function container must have rank 2, 3, 4, or 5.");
1085  INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
1086  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Input field signs container must have rank 2.");
1087  INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.extent(0) != fieldSigns.extent(0), std::invalid_argument,
1088  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Zeroth dimensions (number of integration domains) of the function and field signs containers must agree!");
1089  INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.extent(1) != fieldSigns.extent(1), std::invalid_argument,
1090  ">>> ERROR (FunctionSpaceTools::applyFieldSigns): First dimensions (number of fields) of the function and field signs containers must agree!");
1091 
1092 #endif
1093 
1094  constexpr bool are_accessible =
1095  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1096  typename decltype(inoutFunction)::memory_space>::accessible &&
1097  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1098  typename decltype(fieldSigns)::memory_space>::accessible;
1099  static_assert(are_accessible, "FunctionSpaceTools<DeviceType>::applyFieldSigns(..): input/output views' memory spaces are not compatible with DeviceType");
1100 
1102  <decltype(inoutFunction),decltype(fieldSigns)>;
1103 
1104  const ordinal_type C = inoutFunction.extent(0);
1105  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
1106  Kokkos::parallel_for( policy, FunctorType(inoutFunction, fieldSigns) );
1107  }
1108 
1109  // ------------------------------------------------------------------------------------
1110 
1111  namespace FunctorFunctionSpaceTools {
1112 
1116  template<typename outputPointViewType,
1117  typename inputCoeffViewType,
1118  typename inputFieldViewType>
1119  struct F_evaluate {
1120  outputPointViewType _outputPointVals;
1121  const inputCoeffViewType _inputCoeffs;
1122  const inputFieldViewType _inputFields;
1123 
1124  KOKKOS_INLINE_FUNCTION
1125  F_evaluate( outputPointViewType outputPointVals_,
1126  inputCoeffViewType inputCoeffs_,
1127  inputFieldViewType inputFields_ )
1128  : _outputPointVals(outputPointVals_), _inputCoeffs(inputCoeffs_), _inputFields(inputFields_) {}
1129 
1130  KOKKOS_INLINE_FUNCTION
1131  void operator()(const ordinal_type cl) const {
1132  const ordinal_type nbfs = _inputFields.extent(1);
1133  const ordinal_type npts = _inputFields.extent(2);
1134 
1135  const ordinal_type iend = _inputFields.extent(3);
1136  const ordinal_type jend = _inputFields.extent(4);
1137 
1138  for (ordinal_type bf=0;bf<nbfs;++bf)
1139  for (ordinal_type pt=0;pt<npts;++pt)
1140  for (ordinal_type i=0;i<iend;++i)
1141  for (ordinal_type j=0;j<jend;++j)
1142  _outputPointVals(cl, pt, i, j) += _inputCoeffs(cl, bf) * _inputFields(cl, bf, pt, i, j);
1143  }
1144  };
1145  }
1146 
1147  template<typename DeviceType>
1148  template<typename outputPointValueType, class ...outputPointProperties,
1149  typename inputCoeffValueType, class ...inputCoeffProperties,
1150  typename inputFieldValueType, class ...inputFieldProperties>
1151  void
1153  evaluate( Kokkos::DynRankView<outputPointValueType,outputPointProperties...> outputPointVals,
1154  const Kokkos::DynRankView<inputCoeffValueType, inputCoeffProperties...> inputCoeffs,
1155  const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields ) {
1156 
1157 #ifdef HAVE_INTREPID2_DEBUG
1158  INTREPID2_TEST_FOR_EXCEPTION( inputFields.rank() < 3 || inputFields.rank() > 5, std::invalid_argument,
1159  ">>> ERROR (FunctionSpaceTools::evaluate): Input fields container must have rank 3, 4, or 5.");
1160  INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.rank() != 2, std::invalid_argument,
1161  ">>> ERROR (FunctionSpaceTools::evaluate): Input coefficient container must have rank 2.");
1162  INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.rank() != (inputFields.rank()-1), std::invalid_argument,
1163  ">>> ERROR (FunctionSpaceTools::evaluate): Output values container must have rank one less than the rank of the input fields container.");
1164  INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.extent(0) != inputFields.extent(0), std::invalid_argument,
1165  ">>> ERROR (FunctionSpaceTools::evaluate): Zeroth dimensions (number of cells) of the coefficient and fields input containers must agree!");
1166  INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.extent(1) != inputFields.extent(1), std::invalid_argument,
1167  ">>> ERROR (FunctionSpaceTools::evaluate): First dimensions (number of fields) of the coefficient and fields input containers must agree!");
1168  INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.extent(0) != inputFields.extent(0), std::invalid_argument,
1169  ">>> ERROR (FunctionSpaceTools::evaluate): Zeroth dimensions (number of cells) of the input fields container and the output values container must agree!");
1170  for (size_type i=1;i<outputPointVals.rank();++i)
1171  INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.extent(i) != inputFields.extent(i+1), std::invalid_argument,
1172  ">>> ERROR (FunctionSpaceTools::evaluate): outputPointVals dimension(i) does not match to inputFields dimension(i+1).");
1173 #endif
1174 
1175  constexpr bool are_accessible =
1176  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1177  typename decltype(outputPointVals)::memory_space>::accessible &&
1178  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1179  typename decltype(inputCoeffs)::memory_space>::accessible &&
1180  Kokkos::Impl::MemorySpaceAccess<MemSpaceType,
1181  typename decltype(inputFields)::memory_space>::accessible;
1182  static_assert(are_accessible, "FunctionSpaceTools<DeviceType>::evaluate(..): input/output views' memory spaces are not compatible with DeviceType");
1183 
1184  using FunctorType = FunctorFunctionSpaceTools::F_evaluate
1185  <decltype(outputPointVals),decltype(inputCoeffs),decltype(inputFields)>;
1186 
1187  const ordinal_type C = inputFields.extent(0);
1188  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
1189  Kokkos::parallel_for( policy, FunctorType(outputPointVals, inputCoeffs, inputFields) );
1190  }
1191 
1192  // ------------------------------------------------------------------------------------
1193 
1194 } // end namespace Intrepid2
1195 
1196 #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...
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...
Functor to evaluate functions, see Intrepid2::FunctionSpaceTools for more.
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...