GEOS
H1_Hexahedron_Lagrange1_GaussLegendre2.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 Total, S.A
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
20 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_TRILINEARHEXAHEDRON_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_TRILINEARHEXAHEDRON_HPP_
22 
23 #include "FiniteElementBase.hpp"
24 #include "LagrangeBasis1.hpp"
25 
26 #include <utility>
27 
28 
29 
30 namespace geos
31 {
32 namespace finiteElement
33 {
34 
35 //constexpr static real64 linearBasisAtQuadrature[2] = { 0.5 + 0.5 * 0.5773502691896257645092,
36 // 0.5 - 0.5 * 0.5773502691896257645092 };
37 //__constant__ real64 psiProduct[3] = { 0.5 * linearBasisAtQuadrature[0]*linearBasisAtQuadrature[0],
38 // 0.5 * linearBasisAtQuadrature[0]*linearBasisAtQuadrature[1],
39 // 0.5 * linearBasisAtQuadrature[1]*linearBasisAtQuadrature[1] };
40 
41 //__constant__ real64 psiProduct[3] = { 0.311004233964073108,
42 // 0.083333333333333333,
43 // 0.022329099369260226 };
44 //
45 //__constant__ short dpsi[2] = { -1, 1 };
46 
75 {
76 public:
79 
82 
84  constexpr static localIndex maxSupportPoints = numNodes;
85 
87  constexpr static localIndex numQuadraturePoints = 8;
88 
91 
98  {}
99 
101  virtual localIndex getNumQuadraturePoints() const override
102  {
103  return numQuadraturePoints;
104  }
105 
113  {
114  GEOS_UNUSED_VAR( stack );
115  return numQuadraturePoints;
116  }
117 
119  virtual localIndex getNumSupportPoints() const override
120  {
121  return numNodes;
122  }
123 
125  virtual localIndex getMaxSupportPoints() const override
126  {
127  return maxSupportPoints;
128  }
129 
137  {
138  GEOS_UNUSED_VAR( stack );
139  return numNodes;
140  }
141 
150  static void getSamplingPointCoordInParentSpace( int const linearIndex,
151  real64 (& samplingPointCoord)[3] )
152  {
153  const int i0 = linearIndex % numSamplingPointsPerDirection;
154  const int i1 = ( (linearIndex - i0)/numSamplingPointsPerDirection ) % numSamplingPointsPerDirection;
155  const int i2 = ( (linearIndex - i0)/numSamplingPointsPerDirection - i1 ) / numSamplingPointsPerDirection;
156 
157  constexpr real64 step = 2. / ( numSamplingPointsPerDirection - 1 );
158 
159  samplingPointCoord[0] = -1 + i0 * step;
160  samplingPointCoord[1] = -1 + i1 * step;
161  samplingPointCoord[2] = -1 + i2 * step;
162  }
163 
173  static void calcN( real64 const (&pointCoord)[3],
174  real64 (& N)[numNodes] )
175  {
177  }
178 
187  inline
188  static void calcN( localIndex const q,
189  real64 (& N)[numNodes] )
190  {
191  int qa, qb, qc;
193  real64 const qCoords[3] = { quadratureFactor * LagrangeBasis1::parentSupportCoord( qa ),
194  quadratureFactor * LagrangeBasis1::parentSupportCoord( qb ),
195  quadratureFactor * LagrangeBasis1::parentSupportCoord( qc ) };
196 
197  calcN( qCoords, N );
198  }
199 
209  inline
210  static void calcN( localIndex const q,
211  StackVariables const & stack,
212  real64 ( & N )[numNodes] )
213  {
214  GEOS_UNUSED_VAR( stack );
215  return calcN( q, N );
216  }
217 
227  static void calcFaceBubbleN( real64 const (&pointCoord)[3],
228  real64 (& N)[numFaces] )
229  {
231  }
232 
241  inline
242  static void calcFaceBubbleN( localIndex const q,
243  real64 (& N)[numFaces] )
244  {
245  int qa, qb, qc;
247  real64 const qCoords[3] = { quadratureFactor * LagrangeBasis1::parentSupportCoord( qa ),
248  quadratureFactor * LagrangeBasis1::parentSupportCoord( qb ),
249  quadratureFactor * LagrangeBasis1::parentSupportCoord( qc ) };
250 
251  calcFaceBubbleN( qCoords, N );
252  }
253 
264  static real64 calcGradN( localIndex const q,
265  real64 const (&X)[numNodes][3],
266  real64 ( &gradN )[numNodes][3] );
267 
279  inline
280  static real64 calcGradN( localIndex const q,
281  real64 const (&X)[numNodes][3],
282  StackVariables const & stack,
283  real64 ( &gradN )[numNodes][3] );
284 
296  real64 const (&X)[numNodes][3],
297  real64 ( &gradN )[numFaces][3] );
298 
308  real64 const (&X)[numNodes][3] );
309 
310 
321  real64 const (&X)[numNodes][3],
322  StackVariables const & stack )
323  { GEOS_UNUSED_VAR( stack ); return transformedQuadratureWeight( q, X ); }
324 
325 
335  static real64 invJacobianTransformation( int const q,
336  real64 const (&X)[numNodes][3],
337  real64 ( & J )[3][3] )
338  {
339  int qa, qb, qc;
341  jacobianTransformation( qa, qb, qc, X, J );
342  return LvArray::tensorOps::invert< 3 >( J );
343  }
344 
345 
357  static void symmetricGradient( int const q,
358  real64 const (&invJ)[3][3],
359  real64 const (&var)[numNodes][3],
360  real64 ( &grad )[6] );
361 
362 
363 
380  static void gradient( int const q,
381  real64 const (&invJ)[3][3],
382  real64 const (&var)[numNodes][3],
383  real64 ( &grad )[3][3] );
384 
385 
402  static void plusGradNajAij( int const q,
403  real64 const (&invJ)[3][3],
404  real64 const (&var)[6],
405  real64 ( &R )[numNodes][3] );
406 
407 
408 
419  static void jacobianTransformation( int const qa,
420  int const qb,
421  int const qc,
422  real64 const (&X)[numNodes][3],
423  real64 ( &J )[3][3] );
424 
425 
438  static void
440  int const qb,
441  int const qc,
442  real64 const ( &invJ )[3][3],
443  real64 ( &gradN )[numNodes][3] );
444 
445 
446 private:
448  constexpr static real64 parentLength = LagrangeBasis1::parentSupportCoord( 1 ) - LagrangeBasis1::parentSupportCoord( 0 );
449 
451  constexpr static real64 parentVolume = parentLength*parentLength*parentLength;
452 
454  constexpr static real64 weight = parentVolume / numQuadraturePoints;
455 
459  constexpr static real64 quadratureFactor = 1.0 / 1.732050807568877293528;
460 
461 // constexpr static real64 psiProduct0 = 0.311004233964073108;
462 // constexpr static real64 psiProduct1 = 0.083333333333333333;
463 // constexpr static real64 psiProduct2 = 0.022329099369260226;
464 //
465 // constexpr static short dpsi0 = -1;
466 // constexpr static short dpsi1 = 1;
467 
468 
480  template< typename FUNC, typename ... PARAMS >
482  static void supportLoop( int const qa,
483  int const qb,
484  int const qc,
485  FUNC && func,
486  PARAMS &&... params );
487 
488 };
489 
490 //GEOS_HOST_DEVICE inline real64
491 //psiProductFunc( int const n )
492 //{
493 // // factor = 1./(2 + Sqrt[3]) = 0.267949192431122706
494 //
495 // real64 rval = 0.311004233964073108;
496 // for( int a=0; a<n; ++a )
497 // {
498 // rval *= 0.267949192431122706;
499 // }
500 // return rval;
501 //}
502 
504 
505 template< typename FUNC, typename ... PARAMS >
506 GEOS_HOST_DEVICE inline void
507 H1_Hexahedron_Lagrange1_GaussLegendre2::supportLoop( int const qa,
508  int const qb,
509  int const qc,
510  FUNC && func,
511  PARAMS &&... params )
512 {
513 
515  #define PARENT_GRADIENT_METHOD 2
516 #if PARENT_GRADIENT_METHOD == 1
517  // This option calculates the basis values at the quadrature point for each
518  // linear basis index.
519 
520  real64 const quadratureCoords[3] = { -quadratureFactor + 1.154700538379252 * qa,
521  -quadratureFactor + 1.154700538379252 * qb,
522  -quadratureFactor + 1.154700538379252 * qc };
523 
524  real64 const psi0[2] = { 0.5 - 0.5 * quadratureCoords[0],
525  0.5 + 0.5 * quadratureCoords[0] };
526  real64 const psi1[2] = { 0.5 - 0.5 * quadratureCoords[1],
527  0.5 + 0.5 * quadratureCoords[1] };
528  real64 const psi2[2] = { 0.5 - 0.5 * quadratureCoords[2],
529  0.5 + 0.5 * quadratureCoords[2] };
530  constexpr real64 dpsi[2] = { -0.5, 0.5 };
531 #elif PARENT_GRADIENT_METHOD == 2
532  // This option calculates the product of linear basis prior to use.
533  // The tensor product basis gradient may be expressed as a permutation of the
534  // product between the two possible linear basis gradients. Thus the values
535  // in the basis loop of qaa/qbb/qcc indicate which permutation to choose.
536  // The quantities qaa, qbb, qcc are the difference in index between the
537  // quadrature point and the basis. This is possible because there are 8
538  // basis, and 8 quadrature points, which have correlated indices. So the
539  // values of qaa/qbb/qcc are the "distance" from the quadrature point index
540  // and the support point index.
541  // THIS approach uses about 10 less registers than option 1, with no apparent
542  // cost.
543 
544 // constexpr static real64 linearBasisAtQuadrature[2] = { 0.5 + 0.5 * quadratureFactor,
545 // 0.5 - 0.5 * quadratureFactor };
546 // constexpr static real64 psiProduct[3] = { 0.5 * linearBasisAtQuadrature[0]*linearBasisAtQuadrature[0],
547 // 0.5 * linearBasisAtQuadrature[0]*linearBasisAtQuadrature[1],
548 // 0.5 * linearBasisAtQuadrature[1]*linearBasisAtQuadrature[1] };
549 
551  constexpr static real64 psiProduct[3] = { 0.311004233964073108, 0.083333333333333333, 0.022329099369260226};
552  constexpr static int dpsi[2] = { -1, 1 };
553 
554 // constexpr static real64 psiProduct[3] = { psiProduct0, psiProduct1, psiProduct2 };
555 // constexpr short dpsi[2] = { dpsi0, dpsi1 };
556 #endif
557 
558  // Loop over the linear basis indices in each direction.
559  for( int a=0; a<2; ++a )
560  {
561 #if PARENT_GRADIENT_METHOD == 2
562  int const qaa = ( a^qa ); // abs(a-qa)
563 #endif
564  for( int b=0; b<2; ++b )
565  {
566 #if PARENT_GRADIENT_METHOD == 2
567  int const qbb = ( b^qb );
568 #endif
569  for( int c=0; c<2; ++c )
570  {
571 #if PARENT_GRADIENT_METHOD == 2
572  int const qcc = ( c^qc );
573 #endif
574 
575 #if PARENT_GRADIENT_METHOD == 1
576  real64 const dNdXi[3] = { dpsi[a] * psi1[b] * psi2[c],
577  psi0[a] * dpsi[b] * psi2[c],
578  psi0[a] * psi1[b] * dpsi[c] };
579 #elif PARENT_GRADIENT_METHOD == 2
580 
581 
582 // const real64 dNdXi[3] = { dpsi[a] * psiProductFunc( qbb + qcc ),
583 // dpsi[b] * psiProductFunc( qaa + qcc ),
584 // dpsi[c] * psiProductFunc( qaa + qbb ) };
585 
586  real64 const dNdXi[3] = { dpsi[a] * psiProduct[ qbb + qcc ],
587  dpsi[b] * psiProduct[ qaa + qcc ],
588  dpsi[c] * psiProduct[ qaa + qbb ] };
589 #endif
590  localIndex const nodeIndex = LagrangeBasis1::TensorProduct3D::linearIndex( a, b, c );
591 
592  func( dNdXi, nodeIndex, std::forward< PARAMS >( params )... );
593  }
594  }
595  }
596 }
597 
598 //*************************************************************************************************
600 inline
601 real64
603  real64 const (&X)[numNodes][3],
604  real64 (& gradN)[numNodes][3] )
605 {
606  real64 J[3][3] = {{0}};
607 
608 
609  int qa, qb, qc;
611 
612  jacobianTransformation( qa, qb, qc, X, J );
613 
614  real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
615 
616  applyTransformationToParentGradients( qa, qb, qc, J, gradN );
617 
618  return detJ * weight;
619 }
620 
622 inline
624  calcGradN( localIndex const q,
625  real64 const (&X)[numNodes][3],
626  StackVariables const & GEOS_UNUSED_PARAM( stack ),
627  real64 ( & gradN )[numNodes][3] )
628 {
629  return calcGradN( q, X, gradN );
630 }
631 
633 inline
634 real64
636  real64 const (&X)[numNodes][3],
637  real64 (& gradN)[numFaces][3] )
638 {
639  real64 J[3][3] = {{0}};
640 
641 
642  int qa, qb, qc;
644 
645  jacobianTransformation( qa, qb, qc, X, J );
646 
647  real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
648 
649  real64 dNdXi[numFaces][3] = {{0}};
650 
651  real64 const qCoords[3] = { quadratureFactor * LagrangeBasis1::parentSupportCoord( qa ),
652  quadratureFactor * LagrangeBasis1::parentSupportCoord( qb ),
653  quadratureFactor * LagrangeBasis1::parentSupportCoord( qc ) };
654 
656 
657  for( int fi=0; fi<numFaces; ++fi )
658  {
659  gradN[fi][0] = dNdXi[fi][0] * J[0][0] + dNdXi[fi][1] * J[1][0] + dNdXi[fi][2] * J[2][0];
660  gradN[fi][1] = dNdXi[fi][0] * J[0][1] + dNdXi[fi][1] * J[1][1] + dNdXi[fi][2] * J[2][1];
661  gradN[fi][2] = dNdXi[fi][0] * J[0][2] + dNdXi[fi][1] * J[1][2] + dNdXi[fi][2] * J[2][2];
662  }
663 
664  return detJ * weight;
665 }
666 
667 //*************************************************************************************************
668 #if __GNUC__
669 #pragma GCC diagnostic push
670 #pragma GCC diagnostic ignored "-Wshadow"
671 #endif
672 
674 inline
675 void
677  jacobianTransformation( int const qa,
678  int const qb,
679  int const qc,
680  real64 const (&X)[numNodes][3],
681  real64 ( & J )[3][3] )
682 {
683  supportLoop( qa, qb, qc, [] GEOS_HOST_DEVICE ( real64 const (&dNdXi)[3],
684  int const nodeIndex,
685  real64 const (&X)[numNodes][3],
686  real64 (& J)[3][3] )
687  {
688  real64 const * const GEOS_RESTRICT Xnode = X[nodeIndex];
689  for( int i = 0; i < 3; ++i )
690  {
691  for( int j = 0; j < 3; ++j )
692  {
693  J[i][j] = J[i][j] + dNdXi[ j ] * Xnode[i];
694  }
695  }
696 
697 // J[0][0] = J[0][0] + dNdXi[0] * Xnode[0];
698 // J[0][1] = J[0][1] + dNdXi[1] * Xnode[0];
699 // J[0][2] = J[0][2] + dNdXi[2] * Xnode[0];
700 // J[1][0] = J[1][0] + dNdXi[0] * Xnode[1];
701 // J[1][1] = J[1][1] + dNdXi[1] * Xnode[1];
702 // J[1][2] = J[1][2] + dNdXi[2] * Xnode[1];
703 // J[2][0] = J[2][0] + dNdXi[0] * Xnode[2];
704 // J[2][1] = J[2][1] + dNdXi[1] * Xnode[2];
705 // J[2][2] = J[2][2] + dNdXi[2] * Xnode[2];
706 
707  }, X, J );
708 }
709 
710 
711 //*************************************************************************************************
713 inline
714 void
717  int const qb,
718  int const qc,
719  real64 const ( &invJ )[3][3],
720  real64 (& gradN)[numNodes][3] )
721 {
722  supportLoop( qa, qb, qc, [] GEOS_HOST_DEVICE ( real64 const (&dNdXi)[3],
723  int const nodeIndex,
724  real64 const (&invJ)[3][3],
725  real64 (& gradN)[numNodes][3] )
726  {
727 // for( int i = 0; i < 3; ++i )
728 // {
729 // gradN[nodeIndex][i] = 0.0;
730 // for( int j = 0; j < 3; ++j )
731 // {
732 // gradN[nodeIndex][i] = gradN[nodeIndex][i] + dNdXi[ j ] * invJ[j][i];
733 // }
734 // }
735  // smaller register footprint by manually unrolling the for loops.
736  gradN[nodeIndex][0] = dNdXi[0] * invJ[0][0] + dNdXi[1] * invJ[1][0] + dNdXi[2] * invJ[2][0];
737  gradN[nodeIndex][1] = dNdXi[0] * invJ[0][1] + dNdXi[1] * invJ[1][1] + dNdXi[2] * invJ[2][1];
738  gradN[nodeIndex][2] = dNdXi[0] * invJ[0][2] + dNdXi[1] * invJ[1][2] + dNdXi[2] * invJ[2][2];
739 
740 
741  }, invJ, gradN );
742 }
743 
744 
746 inline
747 real64
750  real64 const (&X)[numNodes][3] )
751 {
752  real64 J[3][3] = {{0}};
753 
754  int qa, qb, qc;
756 
757  jacobianTransformation( qa, qb, qc, X, J );
758 
759  return LvArray::tensorOps::determinant< 3 >( J );
760 }
761 
762 
763 
765 inline
767  real64 const (&invJ)[3][3],
768  real64 const (&var)[numNodes][3],
769  real64 (& grad)[6] )
770 {
771  int qa, qb, qc;
773 
774  supportLoop( qa, qb, qc, [] GEOS_HOST_DEVICE ( real64 const (&dNdXi)[3],
775  int const nodeIndex,
776  real64 const (&invJ)[3][3],
777  real64 const (&var)[numNodes][3],
778  real64 (& grad)[6] )
779  {
780 
781  real64 gradN[3] = {0, 0, 0};
782  for( int i = 0; i < 3; ++i )
783  {
784  for( int j = 0; j < 3; ++j )
785  {
786  gradN[i] = gradN[i] + dNdXi[ j ] * invJ[j][i];
787  }
788  }
789 
790  grad[0] = grad[0] + gradN[0] * var[ nodeIndex ][0];
791  grad[1] = grad[1] + gradN[1] * var[ nodeIndex ][1];
792  grad[2] = grad[2] + gradN[2] * var[ nodeIndex ][2];
793  grad[3] = grad[3] + gradN[2] * var[ nodeIndex ][1] + gradN[1] * var[ nodeIndex ][2];
794  grad[4] = grad[4] + gradN[2] * var[ nodeIndex ][0] + gradN[0] * var[ nodeIndex ][2];
795  grad[5] = grad[5] + gradN[1] * var[ nodeIndex ][0] + gradN[0] * var[ nodeIndex ][1];
796  }, invJ, var, grad );
797 }
798 
800 inline
802  real64 const (&invJ)[3][3],
803  real64 const (&var)[6],
804  real64 (& R)[numNodes][3] )
805 {
806  int qa, qb, qc;
808 
809  supportLoop( qa, qb, qc,
811  ( real64 const (&dNdXi)[3],
812  int const nodeIndex,
813  real64 const (&invJ)[3][3],
814  real64 const (&var)[6],
815  real64 (& R)[numNodes][3] )
816  {
817 
818  real64 gradN[3] = {0, 0, 0};
819  for( int i = 0; i < 3; ++i )
820  {
821  for( int j = 0; j < 3; ++j )
822  {
823  gradN[i] = gradN[i] + dNdXi[ j ] * invJ[j][i];
824  }
825  }
826  R[ nodeIndex ][ 0 ] = R[ nodeIndex ][ 0 ] - var[ 0 ] * gradN[ 0 ] - var[ 5 ] * gradN[ 1 ] - var[ 4 ] * gradN[ 2 ];
827  R[ nodeIndex ][ 1 ] = R[ nodeIndex ][ 1 ] - var[ 5 ] * gradN[ 0 ] - var[ 1 ] * gradN[ 1 ] - var[ 3 ] * gradN[ 2 ];
828  R[ nodeIndex ][ 2 ] = R[ nodeIndex ][ 2 ] - var[ 4 ] * gradN[ 0 ] - var[ 3 ] * gradN[ 1 ] - var[ 2 ] * gradN[ 2 ];
829  }, invJ, var, R );
830 }
831 
832 
833 
835 inline
837  real64 const (&invJ)[3][3],
838  real64 const (&var)[numNodes][3],
839  real64 (& grad)[3][3] )
840 {
841  int qa, qb, qc;
843 
844  supportLoop( qa, qb, qc, [] GEOS_HOST_DEVICE ( real64 const (&dNdXi)[3],
845  int const nodeIndex,
846  real64 const (&invJ)[3][3],
847  real64 const (&var)[numNodes][3],
848  real64 (& grad)[3][3] )
849  {
850  for( int i = 0; i < 3; ++i )
851  {
852  real64 gradN=0.0;;
853  for( int j = 0; j < 3; ++j )
854  {
855  gradN = gradN + dNdXi[ j ] * invJ[j][i];
856  }
857  for( int k = 0; k < 3; ++k )
858  {
859  grad[k][i] = grad[k][i] + gradN * var[ nodeIndex ][k];
860  }
861  }
862  }, invJ, var, grad );
863 }
864 
866 
867 #if __GNUC__
868 #pragma GCC diagnostic pop
869 #endif
870 
871 #undef PARENT_GRADIENT_METHOD
872 }
873 }
874 
875 #endif //FINITE_ELEMENT_SHAPE_KERNEL
#define USING_FINITEELEMENTBASE
Macro to simplify name resolution in derived classes.
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:49
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
#define GEOS_FORCE_INLINE
Marks a function or lambda for inlining.
Definition: GeosxMacros.hpp:51
#define GEOS_RESTRICT
preprocessor variable for the C99 restrict keyword for use with pointers
#define GEOS_UNUSED_PARAM(X)
Mark an unused argument and silence compiler warnings.
Definition: GeosxMacros.hpp:72
Base class for FEM element implementations.
constexpr static int numSamplingPointsPerDirection
Number of sampling points.
constexpr static localIndex maxSupportPoints
The maximum number of support points per element.
static GEOS_HOST_DEVICE void calcN(localIndex const q, real64(&N)[numNodes])
Calculate shape functions values for each support point at a quadrature point.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void getSamplingPointCoordInParentSpace(int const linearIndex, real64(&samplingPointCoord)[3])
Get the Sampling Point Coord In the Parent Space.
static GEOS_HOST_DEVICE real64 calcGradN(localIndex const q, real64 const (&X)[numNodes][3], StackVariables const &stack, real64(&gradN)[numNodes][3])
Calculate the shape functions derivatives wrt the physical coordinates.
virtual GEOS_HOST_DEVICE localIndex getMaxSupportPoints() const override
Get the maximum number of support points for this element.
static GEOS_HOST_DEVICE void calcN(localIndex const q, StackVariables const &stack, real64(&N)[numNodes])
Calculate shape functions values for each support point at a quadrature point.
constexpr static localIndex numNodes
The number of nodes/support points per element.
static GEOS_HOST_DEVICE void applyTransformationToParentGradients(int const qa, int const qb, int const qc, real64 const (&invJ)[3][3], real64(&gradN)[numNodes][3])
Apply a Jacobian transformation matrix from the parent space to the physical space on the parent shap...
virtual GEOS_HOST_DEVICE localIndex getNumSupportPoints() const override
Virtual getter for the number of support points per element.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints(StackVariables const &stack)
Get the number of support points.
static GEOS_HOST_DEVICE void gradient(int const q, real64 const (&invJ)[3][3], real64 const (&var)[numNodes][3], real64(&grad)[3][3])
Calculate the gradient of a vector valued support field at a point using the stored basis function gr...
static GEOS_HOST_DEVICE real64 calcGradN(localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numNodes][3])
Calculate the shape functions derivatives wrt the physical coordinates.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
static GEOS_HOST_DEVICE void jacobianTransformation(int const qa, int const qb, int const qc, real64 const (&X)[numNodes][3], real64(&J)[3][3])
Calculates the isoparametric "Jacobian" transformation matrix/mapping from the parent space to the ph...
static GEOS_HOST_DEVICE void symmetricGradient(int const q, real64 const (&invJ)[3][3], real64 const (&var)[numNodes][3], real64(&grad)[6])
Calculate the symmetric gradient of a vector valued support field at a quadrature point using the sto...
constexpr static int numSamplingPoints
The number of sampling points per element.
static GEOS_HOST_DEVICE real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3], StackVariables const &stack)
Calculate the integration weights for a quadrature point.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcFaceBubbleN(real64 const (&pointCoord)[3], real64(&N)[numFaces])
Calculate face bubble functions values for each face at a given point in the parent space.
virtual GEOS_HOST_DEVICE localIndex getNumQuadraturePoints() const override
Virtual getter for the number of quadrature points per element.
constexpr static localIndex numQuadraturePoints
The number of quadrature points per element.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcN(real64 const (&pointCoord)[3], real64(&N)[numNodes])
Calculate shape functions values for each support point at a given point in the parent space.
static GEOS_HOST_DEVICE real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3])
Calculate the integration weights for a quadrature point.
constexpr static localIndex numFaces
The number of faces/support points for bubble functions per element.
static GEOS_HOST_DEVICE real64 calcGradFaceBubbleN(localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numFaces][3])
Calculate the bubble function derivatives wrt the physical coordinates.
static GEOS_HOST_DEVICE void calcFaceBubbleN(localIndex const q, real64(&N)[numFaces])
Calculate face bubble functions values for each face at a quadrature point.
static GEOS_HOST_DEVICE real64 invJacobianTransformation(int const q, real64 const (&X)[numNodes][3], real64(&J)[3][3])
Calculates the isoparametric "Jacobian" transformation matrix/mapping from the parent space to the ph...
static GEOS_HOST_DEVICE void plusGradNajAij(int const q, real64 const (&invJ)[3][3], real64 const (&var)[6], real64(&R)[numNodes][3])
Inner product of all basis function gradients and a rank-2 symmetric tensor evaluated at a quadrature...
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE real64 parentSupportCoord(const localIndex supportPointIndex)
Calculate the parent coordinates for the xi0 direction, given the linear index of a support point.
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:85
constexpr static localIndex numSupportPoints
The number of support points in the basis.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void value(real64 const (&coords)[3], real64(&N)[numSupportPoints])
The value of the basis function for a support point evaluated at a point along the axes.
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE void multiIndex(const int linearIndex, int &i0, int &i1, int &i2)
Calculate the Cartesian/TensorProduct index given the linear index of a support point.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void gradientFaceBubble(real64 const (&coords)[3], real64(&dNdXi)[numSupportFaces][3])
The value of the bubble basis function derivatives for a support face evaluated at a point along the ...
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE int linearIndex(const int i, const int j, const int k)
Calculates the linear index for support/quadrature points from ijk coordinates.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void valueFaceBubble(real64 const (&coords)[3], real64(&N)[numSupportFaces])
The value of the bubble basis function for a support face evaluated at a point along the axes.
constexpr static localIndex numSupportFaces
The number of support faces in the basis.