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