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 TotalEnergies
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 
81  using StackVariables = typename Base::StackVariables;
82 
84  template< typename SubregionType >
85  using MeshData = typename Base::template MeshData< SubregionType >;
86 
88  using Base::numNodes;
89 
92 
95 
104 
108  {
109  return numQuadraturePoints;
110  }
111 
119  {
120  GEOS_UNUSED_VAR( stack );
121  return numQuadraturePoints;
122  }
123 
127  {
128  return numNodes;
129  }
130 
134  {
135  return maxSupportPoints;
136  }
137 
145  {
146  GEOS_UNUSED_VAR( stack );
147  return numNodes;
148  }
149 
158  static void getSamplingPointCoordInParentSpace( int const linearIndex,
159  real64 (& samplingPointCoord)[3] )
160  {
161  const int i0 = linearIndex % numSamplingPointsPerDirection;
162  const int i1 = ( (linearIndex - i0)/numSamplingPointsPerDirection ) % numSamplingPointsPerDirection;
163  const int i2 = ( (linearIndex - i0)/numSamplingPointsPerDirection - i1 ) / numSamplingPointsPerDirection;
164 
165  constexpr real64 step = 2. / ( numSamplingPointsPerDirection - 1 );
166 
167  samplingPointCoord[0] = -1 + i0 * step;
168  samplingPointCoord[1] = -1 + i1 * step;
169  samplingPointCoord[2] = -1 + i2 * step;
170  }
171 
181  static void calcN( real64 const (&pointCoord)[3],
182  real64 (& N)[numNodes] )
183  {
185  }
186 
195  inline
196  static void calcN( localIndex const q,
197  real64 (& N)[numNodes] )
198  {
199  int qa, qb, qc;
201  real64 const qCoords[3] = { quadratureFactor * LagrangeBasis1::parentSupportCoord( qa ),
202  quadratureFactor * LagrangeBasis1::parentSupportCoord( qb ),
203  quadratureFactor * LagrangeBasis1::parentSupportCoord( qc ) };
204 
205  calcN( qCoords, N );
206  }
207 
217  inline
218  static void calcN( localIndex const q,
219  StackVariables const & stack,
220  real64 ( & N )[numNodes] )
221  {
222  GEOS_UNUSED_VAR( stack );
223  return calcN( q, N );
224  }
225 
235  static void calcFaceBubbleN( real64 const (&pointCoord)[3],
236  real64 (& N)[numFaces] )
237  {
239  }
240 
249  inline
250  static void calcFaceBubbleN( localIndex const q,
251  real64 (& N)[numFaces] )
252  {
253  int qa, qb, qc;
255  real64 const qCoords[3] = { quadratureFactor * LagrangeBasis1::parentSupportCoord( qa ),
256  quadratureFactor * LagrangeBasis1::parentSupportCoord( qb ),
257  quadratureFactor * LagrangeBasis1::parentSupportCoord( qc ) };
258 
259  calcFaceBubbleN( qCoords, N );
260  }
261 
271  static
273  real64 const (&X)[numNodes][3],
274  real64 ( &J )[3][3] );
275 
286  static
288  real64 const (&X)[numNodes][3],
289  StackVariables const & stack,
290  real64 ( &J )[3][3] );
291 
302  static real64 calcGradN( localIndex const q,
303  real64 const (&X)[numNodes][3],
304  real64 ( &gradN )[numNodes][3] );
305 
317  inline
318  static real64 calcGradN( localIndex const q,
319  real64 const (&X)[numNodes][3],
320  StackVariables const & stack,
321  real64 ( &gradN )[numNodes][3] );
322 
334  real64 const (&X)[numNodes][3],
335  real64 ( &gradN )[numFaces][3] );
336 
346  real64 const (&X)[numNodes][3] );
347 
348 
359  real64 const (&X)[numNodes][3],
360  StackVariables const & stack )
361  { GEOS_UNUSED_VAR( stack ); return transformedQuadratureWeight( q, X ); }
362 
363 
373  static real64 invJacobianTransformation( int const q,
374  real64 const (&X)[numNodes][3],
375  real64 ( & J )[3][3] )
376  {
377  int qa, qb, qc;
379  jacobianTransformation( qa, qb, qc, X, J );
380  return LvArray::tensorOps::invert< 3 >( J );
381  }
382 
383 
395  static void symmetricGradient( int const q,
396  real64 const (&invJ)[3][3],
397  real64 const (&var)[numNodes][3],
398  real64 ( &grad )[6] );
399 
400 
401 
418  static void gradient( int const q,
419  real64 const (&invJ)[3][3],
420  real64 const (&var)[numNodes][3],
421  real64 ( &grad )[3][3] );
422 
423 
440  static void plusGradNajAij( int const q,
441  real64 const (&invJ)[3][3],
442  real64 const (&var)[6],
443  real64 ( &R )[numNodes][3] );
444 
445 
446 
457  static void jacobianTransformation( int const qa,
458  int const qb,
459  int const qc,
460  real64 const (&X)[numNodes][3],
461  real64 ( &J )[3][3] );
462 
463 
476  static void
478  int const qb,
479  int const qc,
480  real64 const ( &invJ )[3][3],
481  real64 ( &gradN )[numNodes][3] );
482 
483 
484 private:
486  constexpr static real64 parentLength = LagrangeBasis1::parentSupportCoord( 1 ) - LagrangeBasis1::parentSupportCoord( 0 );
487 
489  constexpr static real64 parentVolume = parentLength*parentLength*parentLength;
490 
492  constexpr static real64 weight = parentVolume / numQuadraturePoints;
493 
497  constexpr static real64 quadratureFactor = 1.0 / 1.732050807568877293528;
498 
499 // constexpr static real64 psiProduct0 = 0.311004233964073108;
500 // constexpr static real64 psiProduct1 = 0.083333333333333333;
501 // constexpr static real64 psiProduct2 = 0.022329099369260226;
502 //
503 // constexpr static short dpsi0 = -1;
504 // constexpr static short dpsi1 = 1;
505 
506 
518  template< typename FUNC, typename ... PARAMS >
520  static void supportLoop( int const qa,
521  int const qb,
522  int const qc,
523  FUNC && func,
524  PARAMS &&... params );
525 
526 };
527 
528 //GEOS_HOST_DEVICE inline real64
529 //psiProductFunc( int const n )
530 //{
531 // // factor = 1./(2 + Sqrt[3]) = 0.267949192431122706
532 //
533 // real64 rval = 0.311004233964073108;
534 // for( int a=0; a<n; ++a )
535 // {
536 // rval *= 0.267949192431122706;
537 // }
538 // return rval;
539 //}
540 
542 
543 template< typename FUNC, typename ... PARAMS >
544 GEOS_HOST_DEVICE inline void
545 H1_Hexahedron_Lagrange1_GaussLegendre2_impl::supportLoop( int const qa,
546  int const qb,
547  int const qc,
548  FUNC && func,
549  PARAMS &&... params )
550 {
551 
553  #define PARENT_GRADIENT_METHOD 2
554 #if PARENT_GRADIENT_METHOD == 1
555  // This option calculates the basis values at the quadrature point for each
556  // linear basis index.
557 
558  real64 const quadratureCoords[3] = { -quadratureFactor + 1.154700538379252 * qa,
559  -quadratureFactor + 1.154700538379252 * qb,
560  -quadratureFactor + 1.154700538379252 * qc };
561 
562  real64 const psi0[2] = { 0.5 - 0.5 * quadratureCoords[0],
563  0.5 + 0.5 * quadratureCoords[0] };
564  real64 const psi1[2] = { 0.5 - 0.5 * quadratureCoords[1],
565  0.5 + 0.5 * quadratureCoords[1] };
566  real64 const psi2[2] = { 0.5 - 0.5 * quadratureCoords[2],
567  0.5 + 0.5 * quadratureCoords[2] };
568  constexpr real64 dpsi[2] = { -0.5, 0.5 };
569 #elif PARENT_GRADIENT_METHOD == 2
570  // This option calculates the product of linear basis prior to use.
571  // The tensor product basis gradient may be expressed as a permutation of the
572  // product between the two possible linear basis gradients. Thus the values
573  // in the basis loop of qaa/qbb/qcc indicate which permutation to choose.
574  // The quantities qaa, qbb, qcc are the difference in index between the
575  // quadrature point and the basis. This is possible because there are 8
576  // basis, and 8 quadrature points, which have correlated indices. So the
577  // values of qaa/qbb/qcc are the "distance" from the quadrature point index
578  // and the support point index.
579  // THIS approach uses about 10 less registers than option 1, with no apparent
580  // cost.
581 
582 // constexpr static real64 linearBasisAtQuadrature[2] = { 0.5 + 0.5 * quadratureFactor,
583 // 0.5 - 0.5 * quadratureFactor };
584 // constexpr static real64 psiProduct[3] = { 0.5 * linearBasisAtQuadrature[0]*linearBasisAtQuadrature[0],
585 // 0.5 * linearBasisAtQuadrature[0]*linearBasisAtQuadrature[1],
586 // 0.5 * linearBasisAtQuadrature[1]*linearBasisAtQuadrature[1] };
587 
589  constexpr static real64 psiProduct[3] = { 0.311004233964073108, 0.083333333333333333, 0.022329099369260226};
590  constexpr static int dpsi[2] = { -1, 1 };
591 
592 // constexpr static real64 psiProduct[3] = { psiProduct0, psiProduct1, psiProduct2 };
593 // constexpr short dpsi[2] = { dpsi0, dpsi1 };
594 #endif
595 
596  // Loop over the linear basis indices in each direction.
597  for( int a=0; a<2; ++a )
598  {
599 #if PARENT_GRADIENT_METHOD == 2
600  int const qaa = ( a^qa ); // abs(a-qa)
601 #endif
602  for( int b=0; b<2; ++b )
603  {
604 #if PARENT_GRADIENT_METHOD == 2
605  int const qbb = ( b^qb );
606 #endif
607  for( int c=0; c<2; ++c )
608  {
609 #if PARENT_GRADIENT_METHOD == 2
610  int const qcc = ( c^qc );
611 #endif
612 
613 #if PARENT_GRADIENT_METHOD == 1
614  real64 const dNdXi[3] = { dpsi[a] * psi1[b] * psi2[c],
615  psi0[a] * dpsi[b] * psi2[c],
616  psi0[a] * psi1[b] * dpsi[c] };
617 #elif PARENT_GRADIENT_METHOD == 2
618 
619 
620 // const real64 dNdXi[3] = { dpsi[a] * psiProductFunc( qbb + qcc ),
621 // dpsi[b] * psiProductFunc( qaa + qcc ),
622 // dpsi[c] * psiProductFunc( qaa + qbb ) };
623 
624  real64 const dNdXi[3] = { dpsi[a] * psiProduct[ qbb + qcc ],
625  dpsi[b] * psiProduct[ qaa + qcc ],
626  dpsi[c] * psiProduct[ qaa + qbb ] };
627 #endif
628  localIndex const nodeIndex = LagrangeBasis1::TensorProduct3D::linearIndex( a, b, c );
629 
630  func( dNdXi, nodeIndex, std::forward< PARAMS >( params )... );
631  }
632  }
633  }
634 }
635 
636 //*************************************************************************************************
639 real64
641  real64 const (&X)[numNodes][3],
642  real64 (& J)[3][3] )
643 {
644  for( int i = 0; i < 3; ++i )
645  {
646  for( int j = 0; j < 3; ++j )
647  {
648  J[i][j] = 0;
649  }
650  }
651  int qa, qb, qc;
653  jacobianTransformation( qa, qb, qc, X, J );
654  real64 const detJ = LvArray::tensorOps::determinant< 3 >( J );
655  return detJ*weight;
656 }
657 
660 real64
662  real64 const (&X)[numNodes][3],
663  StackVariables const & GEOS_UNUSED_PARAM( stack ),
664  real64 (& J)[3][3] )
665 {
666  return calcJacobian( q, X, J );
667 }
668 
670 inline
671 real64
673  real64 const (&X)[numNodes][3],
674  real64 (& gradN)[numNodes][3] )
675 {
676  real64 J[3][3] = {{0}};
677 
678  int qa, qb, qc;
680 
681  jacobianTransformation( qa, qb, qc, X, J );
682 
683  real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
684 
685  applyTransformationToParentGradients( qa, qb, qc, J, gradN );
686 
687  return detJ * weight;
688 }
689 
691 inline
693  calcGradN( localIndex const q,
694  real64 const (&X)[numNodes][3],
695  StackVariables const & GEOS_UNUSED_PARAM( stack ),
696  real64 ( & gradN )[numNodes][3] )
697 {
698  return calcGradN( q, X, gradN );
699 }
700 
702 inline
703 real64
705  real64 const (&X)[numNodes][3],
706  real64 (& gradN)[numFaces][3] )
707 {
708  real64 J[3][3] = {{0}};
709 
710 
711  int qa, qb, qc;
713 
714  jacobianTransformation( qa, qb, qc, X, J );
715 
716  real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
717 
718  real64 dNdXi[numFaces][3] = {{0}};
719 
720  real64 const qCoords[3] = { quadratureFactor * LagrangeBasis1::parentSupportCoord( qa ),
721  quadratureFactor * LagrangeBasis1::parentSupportCoord( qb ),
722  quadratureFactor * LagrangeBasis1::parentSupportCoord( qc ) };
723 
725 
726  for( int fi=0; fi<numFaces; ++fi )
727  {
728  gradN[fi][0] = dNdXi[fi][0] * J[0][0] + dNdXi[fi][1] * J[1][0] + dNdXi[fi][2] * J[2][0];
729  gradN[fi][1] = dNdXi[fi][0] * J[0][1] + dNdXi[fi][1] * J[1][1] + dNdXi[fi][2] * J[2][1];
730  gradN[fi][2] = dNdXi[fi][0] * J[0][2] + dNdXi[fi][1] * J[1][2] + dNdXi[fi][2] * J[2][2];
731  }
732 
733  return detJ * weight;
734 }
735 
736 //*************************************************************************************************
737 #if __GNUC__
738 #pragma GCC diagnostic push
739 #pragma GCC diagnostic ignored "-Wshadow"
740 #endif
741 
743 inline
744 void
746  jacobianTransformation( int const qa,
747  int const qb,
748  int const qc,
749  real64 const (&X)[numNodes][3],
750  real64 ( & J )[3][3] )
751 {
752  supportLoop( qa, qb, qc, [] GEOS_HOST_DEVICE ( real64 const (&dNdXi)[3],
753  int const nodeIndex,
754  real64 const (&X)[numNodes][3],
755  real64 (& J)[3][3] )
756  {
757  real64 const * const GEOS_RESTRICT Xnode = X[nodeIndex];
758  for( int i = 0; i < 3; ++i )
759  {
760  for( int j = 0; j < 3; ++j )
761  {
762  J[i][j] = J[i][j] + dNdXi[ j ] * Xnode[i];
763  }
764  }
765 
766 // J[0][0] = J[0][0] + dNdXi[0] * Xnode[0];
767 // J[0][1] = J[0][1] + dNdXi[1] * Xnode[0];
768 // J[0][2] = J[0][2] + dNdXi[2] * Xnode[0];
769 // J[1][0] = J[1][0] + dNdXi[0] * Xnode[1];
770 // J[1][1] = J[1][1] + dNdXi[1] * Xnode[1];
771 // J[1][2] = J[1][2] + dNdXi[2] * Xnode[1];
772 // J[2][0] = J[2][0] + dNdXi[0] * Xnode[2];
773 // J[2][1] = J[2][1] + dNdXi[1] * Xnode[2];
774 // J[2][2] = J[2][2] + dNdXi[2] * Xnode[2];
775 
776  }, X, J );
777 }
778 
779 
780 //*************************************************************************************************
782 inline
783 void
786  int const qb,
787  int const qc,
788  real64 const ( &invJ )[3][3],
789  real64 (& gradN)[numNodes][3] )
790 {
791  supportLoop( qa, qb, qc, [] GEOS_HOST_DEVICE ( real64 const (&dNdXi)[3],
792  int const nodeIndex,
793  real64 const (&invJ)[3][3],
794  real64 (& gradN)[numNodes][3] )
795  {
796 // for( int i = 0; i < 3; ++i )
797 // {
798 // gradN[nodeIndex][i] = 0.0;
799 // for( int j = 0; j < 3; ++j )
800 // {
801 // gradN[nodeIndex][i] = gradN[nodeIndex][i] + dNdXi[ j ] * invJ[j][i];
802 // }
803 // }
804  // smaller register footprint by manually unrolling the for loops.
805  gradN[nodeIndex][0] = dNdXi[0] * invJ[0][0] + dNdXi[1] * invJ[1][0] + dNdXi[2] * invJ[2][0];
806  gradN[nodeIndex][1] = dNdXi[0] * invJ[0][1] + dNdXi[1] * invJ[1][1] + dNdXi[2] * invJ[2][1];
807  gradN[nodeIndex][2] = dNdXi[0] * invJ[0][2] + dNdXi[1] * invJ[1][2] + dNdXi[2] * invJ[2][2];
808 
809 
810  }, invJ, gradN );
811 }
812 
813 
815 inline
816 real64
819  real64 const (&X)[numNodes][3] )
820 {
821  real64 J[3][3] = {{0}};
822 
823  int qa, qb, qc;
825 
826  jacobianTransformation( qa, qb, qc, X, J );
827 
828  return LvArray::tensorOps::determinant< 3 >( J );
829 }
830 
831 
832 
834 inline
836  real64 const (&invJ)[3][3],
837  real64 const (&var)[numNodes][3],
838  real64 (& grad)[6] )
839 {
840  int qa, qb, qc;
842 
843  supportLoop( qa, qb, qc, [] GEOS_HOST_DEVICE ( real64 const (&dNdXi)[3],
844  int const nodeIndex,
845  real64 const (&invJ)[3][3],
846  real64 const (&var)[numNodes][3],
847  real64 (& grad)[6] )
848  {
849 
850  real64 gradN[3] = {0, 0, 0};
851  for( int i = 0; i < 3; ++i )
852  {
853  for( int j = 0; j < 3; ++j )
854  {
855  gradN[i] = gradN[i] + dNdXi[ j ] * invJ[j][i];
856  }
857  }
858 
859  grad[0] = grad[0] + gradN[0] * var[ nodeIndex ][0];
860  grad[1] = grad[1] + gradN[1] * var[ nodeIndex ][1];
861  grad[2] = grad[2] + gradN[2] * var[ nodeIndex ][2];
862  grad[3] = grad[3] + gradN[2] * var[ nodeIndex ][1] + gradN[1] * var[ nodeIndex ][2];
863  grad[4] = grad[4] + gradN[2] * var[ nodeIndex ][0] + gradN[0] * var[ nodeIndex ][2];
864  grad[5] = grad[5] + gradN[1] * var[ nodeIndex ][0] + gradN[0] * var[ nodeIndex ][1];
865  }, invJ, var, grad );
866 }
867 
869 inline
871  real64 const (&invJ)[3][3],
872  real64 const (&var)[6],
873  real64 (& R)[numNodes][3] )
874 {
875  int qa, qb, qc;
877 
878  supportLoop( qa, qb, qc,
880  ( real64 const (&dNdXi)[3],
881  int const nodeIndex,
882  real64 const (&invJ)[3][3],
883  real64 const (&var)[6],
884  real64 (& R)[numNodes][3] )
885  {
886 
887  real64 gradN[3] = {0, 0, 0};
888  for( int i = 0; i < 3; ++i )
889  {
890  for( int j = 0; j < 3; ++j )
891  {
892  gradN[i] = gradN[i] + dNdXi[ j ] * invJ[j][i];
893  }
894  }
895  R[ nodeIndex ][ 0 ] = R[ nodeIndex ][ 0 ] - var[ 0 ] * gradN[ 0 ] - var[ 5 ] * gradN[ 1 ] - var[ 4 ] * gradN[ 2 ];
896  R[ nodeIndex ][ 1 ] = R[ nodeIndex ][ 1 ] - var[ 5 ] * gradN[ 0 ] - var[ 1 ] * gradN[ 1 ] - var[ 3 ] * gradN[ 2 ];
897  R[ nodeIndex ][ 2 ] = R[ nodeIndex ][ 2 ] - var[ 4 ] * gradN[ 0 ] - var[ 3 ] * gradN[ 1 ] - var[ 2 ] * gradN[ 2 ];
898  }, invJ, var, R );
899 }
900 
901 
902 
904 inline
906  real64 const (&invJ)[3][3],
907  real64 const (&var)[numNodes][3],
908  real64 (& grad)[3][3] )
909 {
910  int qa, qb, qc;
912 
913  supportLoop( qa, qb, qc, [] GEOS_HOST_DEVICE ( real64 const (&dNdXi)[3],
914  int const nodeIndex,
915  real64 const (&invJ)[3][3],
916  real64 const (&var)[numNodes][3],
917  real64 (& grad)[3][3] )
918  {
919  for( int i = 0; i < 3; ++i )
920  {
921  real64 gradN=0.0;;
922  for( int j = 0; j < 3; ++j )
923  {
924  gradN = gradN + dNdXi[ j ] * invJ[j][i];
925  }
926  for( int k = 0; k < 3; ++k )
927  {
928  grad[k][i] = grad[k][i] + gradN * var[ nodeIndex ][k];
929  }
930  }
931  }, invJ, var, grad );
932 }
933 
935 
936 #if __GNUC__
937 #pragma GCC diagnostic pop
938 #endif
939 
942 {
943 public:
944 
947 
950 
955  { }
956 
958  virtual ~H1_Hexahedron_Lagrange1_GaussLegendre2() override final = default;
959 
965  {
966  return static_cast< H1_Hexahedron_Lagrange1_GaussLegendre2_impl * >(this);
967  }
968 
974  {
975  return static_cast< H1_Hexahedron_Lagrange1_GaussLegendre2_impl const * >(this);
976  }
977 
978 };
979 
980 #undef PARENT_GRADIENT_METHOD
981 }
982 }
983 
984 #endif //FINITE_ELEMENT_SHAPE_KERNEL
#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
Device-compatible (non virtual) Base class for all finite element formulations.
constexpr static localIndex numQuadraturePoints
The number of quadrature points per element.
constexpr static int numSamplingPointsPerDirection
Number of sampling points.
constexpr static localIndex maxSupportPoints
The maximum number of support points per element.
constexpr static localIndex numNodes
The number of nodes per element.
constexpr static localIndex numFaces
The number of faces per element.
GEOS_HOST_DEVICE FiniteElementBase_impl & operator=(FiniteElementBase_impl const &)=default
Default copy assignment operator.
Base class for FEM element implementations.
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...
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.
DO_NOT_DOCUMENT static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints()
Get the number of quadrature points.
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.
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.
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 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 real64 calcJacobian(localIndex const q, real64 const (&X)[numNodes][3], StackVariables const &stack, real64(&J)[3][3])
Calculate the Jacobian matrix at a quadrature point.
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.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void getSamplingPointCoordInParentSpace(int const linearIndex, real64(&samplingPointCoord)[3])
Get the Sampling Point Coord In the Parent Space.
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 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...
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 void calcN(localIndex const q, real64(&N)[numNodes])
Calculate shape functions values for each support point at a quadrature point.
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 localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
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...
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 calcJacobian(localIndex const q, real64 const (&X)[numNodes][3], real64(&J)[3][3])
Calculate the Jacobian matrix at a quadrature point.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints()
Get the number of support points.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints(StackVariables const &stack)
Get the number of support points.
typename Base::template MeshData< SubregionType > MeshData
Mesh data structure for the element.
static GEOS_HOST_DEVICE localIndex getMaxSupportPoints()
Get the maximum number of support points.
static GEOS_HOST_DEVICE real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3])
Calculate the integration weights for 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 real64 calcGradFaceBubbleN(localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numFaces][3])
Calculate the bubble function derivatives wrt the physical coordinates.
H1_Hexahedron_Lagrange1_GaussLegendre2_impl const * getImpl() const
Get the device-compatible implementation type.
H1_Hexahedron_Lagrange1_GaussLegendre2_impl * getImpl()
Get the device-compatible implementation type.
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:98
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:84
Kernel variables (dof numbers, jacobian and residual) located on the stack.
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.