GEOS
H1_Tetrahedron_Lagrange1_Gauss.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_H1TETRAHEDRONLAGRANGE1GAUSS_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1TETRAHEDRONLAGRANGE1GAUSS_HPP_
22 
23 #include "FiniteElementBase.hpp"
24 #include "LagrangeBasis1.hpp"
25 
26 
27 namespace geos
28 {
29 namespace finiteElement
30 {
31 
49 template< typename NUM_Q_POINTS >
50 class H1_Tetrahedron_Lagrange1_Gauss_impl : public FiniteElementBase_impl< 4, 4, NUM_Q_POINTS::value >
51 {
52 public:
55 
57  using StackVariables = typename Base::StackVariables;
58 
60  template< typename SubregionType >
61  using MeshData = typename Base::template MeshData< SubregionType >;
62 
64  using Base::numNodes;
65 
68 
71 
74 
76  constexpr static localIndex numFaces = Base::numFaces;
77 
79  constexpr static int numSamplingPointsPerDirection = 10;
80 
83 
85  static_assert( ( NUM_Q_POINTS::value == 1 ||
86  NUM_Q_POINTS::value == 5 ||
87  NUM_Q_POINTS::value == 14 ),
88  "NUM_Q_POINTS::value must be 1, 5, or 14!" );
89 
92 
101 
105  {
106  return numQuadraturePoints;
107  }
108 
115  static localIndex
117  {
118  GEOS_UNUSED_VAR( stack );
119  return numQuadraturePoints;
120  }
121 
125  {
126  return numNodes;
127  }
128 
132  {
133  return maxSupportPoints;
134  }
135 
143  {
144  GEOS_UNUSED_VAR( stack );
145  return numNodes;
146  }
147 
156  static void getSamplingPointCoordInParentSpace( int const & linearIndex,
157  real64 (& samplingPointCoord)[3] )
158  {
159  int const i0 = linearIndex % numSamplingPointsPerDirection;
160  int const i1 = ( (linearIndex - i0)/numSamplingPointsPerDirection ) % numSamplingPointsPerDirection;
161  int const i2 = ( (linearIndex - i0)/numSamplingPointsPerDirection - i1 ) / numSamplingPointsPerDirection;
162 
163  real64 const step = 1. / ( numSamplingPointsPerDirection - 1 );
164 
165  real64 const r = i0 * step;
166  real64 const s = i1 * step;
167  real64 const t = i2 * step;
168  if( (r+s) <= t )
169  {
170  samplingPointCoord[0] = r;
171  samplingPointCoord[1] = s;
172  samplingPointCoord[2] = t;
173  }
174  else
175  {
176  // if outside of the triangle need to reproject it. Points will be doubled though.
177  samplingPointCoord[0] = 1 - t - r;
178  samplingPointCoord[1] = 1 - t - s;
179  samplingPointCoord[2] = t;
180  }
181  }
182 
192  static void calcN( real64 const (&pointCoord)[3],
193  real64 ( &N )[numNodes] );
194 
203  static void calcN( localIndex const q,
204  real64 ( &N )[numNodes] );
205 
215  inline
216  static void calcN( localIndex const q,
217  StackVariables const & stack,
218  real64 ( &N )[numNodes] );
219 
229  static void calcFaceBubbleN( real64 const (&pointCoord)[3],
230  real64 (& N)[numFaces] )
231  {
232 
233  real64 const r = pointCoord[0];
234  real64 const s = pointCoord[1];
235  real64 const t = pointCoord[2];
236 
237  N[0] = (1 - r - s - t) * r * t;
238  N[1] = (1 - r - s - t) * r * s;
239  N[2] = (1 - r - s - t) * t * s;
240  N[3] = r * s * t;
241  }
242 
251  inline
252  static void calcFaceBubbleN( localIndex const q,
253  real64 (& N)[numFaces] )
254  {
255 
256  real64 const pointCoord[3] = {quadratureParentCoords0( q ),
257  quadratureParentCoords1( q ),
258  quadratureParentCoords2( q )};
259 
260  calcFaceBubbleN( pointCoord, N );
261  }
262 
272  static
274  real64 const (&X)[numNodes][3],
275  real64 ( &J )[3][3] );
276 
287  static
289  real64 const (&X)[numNodes][3],
290  StackVariables const &stack,
291  real64 ( &J )[3][3] );
292 
303  static real64 calcGradN( localIndex const q,
304  real64 const (&X)[numNodes][3],
305  real64 ( &gradN )[numNodes][3] );
306 
318  inline
319  static real64 calcGradN( localIndex const q,
320  real64 const (&X)[numNodes][3],
321  StackVariables const & stack,
322  real64 ( &gradN )[numNodes][3] );
323 
335  real64 const (&X)[numNodes][3],
336  real64 ( &gradN )[numFaces][3] );
337 
347  real64 const (&X)[numNodes][3] );
348 
359  real64 const (&X)[numNodes][3],
360  StackVariables const & stack )
361  { GEOS_UNUSED_VAR( stack ); return transformedQuadratureWeight( q, X ); }
362 
372  static real64 invJacobianTransformation( int const q,
373  real64 const (&X)[numNodes][3],
374  real64 ( & J )[3][3] )
375  {
376  GEOS_UNUSED_VAR( q, X );
377 // jacobianTransformation( q, X, J );
378  return LvArray::tensorOps::invert< 3 >( J );
379  }
380 private:
382  constexpr static real64 parentVolume = 1.0 / 6.0;
383 
385  constexpr static real64 weight = parentVolume;
386  //constexpr static real64 weight = parentVolume / numQuadraturePoints;
387 
394  inline
395  constexpr static real64 quadratureWeight( localIndex const q )
396  {
397  if constexpr (numQuadraturePoints == 1)
398  {
399  real64 const w[numQuadraturePoints] = { 1.0 };
400  return w[q];
401  }
402  else if constexpr (numQuadraturePoints == 5)
403  {
404  real64 const w[numQuadraturePoints] = {-4.0/5.0, 9.0/20.0, 9.0/20.0, 9.0/20.0, 9.0/20.0 };
405  return w[q];
406  }
407  else if constexpr (numQuadraturePoints == 14)
408  {
409  real64 const w[numQuadraturePoints] = { 0.073493043116361949544,
410  0.073493043116361949544,
411  0.073493043116361949544,
412  0.073493043116361949544,
413  0.11268792571801585080,
414  0.11268792571801585080,
415  0.11268792571801585080,
416  0.11268792571801585080,
417  0.042546020777081466438,
418  0.042546020777081466438,
419  0.042546020777081466438,
420  0.042546020777081466438,
421  0.042546020777081466438,
422  0.042546020777081466438 };
423  return w[q];
424  }
425 
426  }
427 
435  inline
436  constexpr static real64 quadratureParentCoords0( localIndex const q )
437  {
438 
439  if constexpr (numQuadraturePoints == 1)
440  {
441  real64 const qCoords[numQuadraturePoints] = { 1.0/4.0 };
442  return qCoords[q];
443  }
444  else if constexpr (numQuadraturePoints == 5)
445  {
446  real64 const qCoords[numQuadraturePoints] = { 1.0/4.0, 1.0/2.0, 1.0/6.0, 1.0/6.0, 1.0/6.0 };
447  return qCoords[q];
448  }
449  else if constexpr (numQuadraturePoints == 14)
450  {
451  real64 const qCoords[numQuadraturePoints] = { 0.72179424906732632079,
452  0.092735250310891226402,
453  0.092735250310891226402,
454  0.092735250310891226402,
455  0.067342242210098170608,
456  0.31088591926330060980,
457  0.31088591926330060980,
458  0.31088591926330060980,
459  0.045503704125649649492,
460  0.045503704125649649492,
461  0.045503704125649649492,
462  0.45449629587435035051,
463  0.45449629587435035051,
464  0.45449629587435035051 };
465  return qCoords[q];
466  }
467 
468  }
469 
477  inline
478  constexpr static real64 quadratureParentCoords1( localIndex const q )
479  {
480 
481  if constexpr (numQuadraturePoints == 1)
482  {
483  real64 const qCoords[numQuadraturePoints] = { 1.0/4.0 };
484  return qCoords[q];
485  }
486  else if constexpr (numQuadraturePoints == 5)
487  {
488  real64 const qCoords[numQuadraturePoints] = { 1.0/4.0, 1.0/6.0, 1.0/2.0, 1.0/6.0, 1.0/6.0 };
489  return qCoords[q];
490  }
491  else if constexpr (numQuadraturePoints == 14)
492  {
493  real64 const qCoords[numQuadraturePoints] = { 0.092735250310891226402,
494  0.72179424906732632079,
495  0.092735250310891226402,
496  0.092735250310891226402,
497  0.31088591926330060980,
498  0.067342242210098170608,
499  0.31088591926330060980,
500  0.31088591926330060980,
501  0.045503704125649649492,
502  0.45449629587435035051,
503  0.45449629587435035051,
504  0.045503704125649649492,
505  0.045503704125649649492,
506  0.45449629587435035051 };
507  return qCoords[q];
508  }
509 
510  }
511 
519  inline
520  constexpr static real64 quadratureParentCoords2( localIndex const q )
521  {
522 
523  if constexpr (numQuadraturePoints == 1)
524  {
525  real64 const qCoords[numQuadraturePoints] = { 1.0/4.0 };
526  return qCoords[q];
527  }
528  else if constexpr (numQuadraturePoints == 5)
529  {
530  real64 const qCoords[numQuadraturePoints] = { 1.0/4.0, 1.0/6.0, 1.0/6.0, 1.0/2.0, 1.0/6.0 };
531  return qCoords[q];
532  }
533  else if constexpr (numQuadraturePoints == 14)
534  {
535  real64 const qCoords[numQuadraturePoints] = { 0.092735250310891226402,
536  0.092735250310891226402,
537  0.72179424906732632079,
538  0.092735250310891226402,
539  0.31088591926330060980,
540  0.31088591926330060980,
541  0.067342242210098170608,
542  0.31088591926330060980,
543  0.45449629587435035051,
544  0.045503704125649649492,
545  0.45449629587435035051,
546  0.045503704125649649492,
547  0.45449629587435035051,
548  0.045503704125649649492 };
549 
550  return qCoords[q];
551  }
552 
553  }
554 
562  static real64 determinantJacobianTransformation( real64 const (&X)[numNodes][3] );
563 
564 };
565 
567 
568 template< typename NUM_Q_POINTS >
570 inline
571 real64
572 H1_Tetrahedron_Lagrange1_Gauss_impl< NUM_Q_POINTS >::
573 determinantJacobianTransformation( real64 const (&X)[numNodes][3] )
574 {
575  return ( X[1][0] - X[0][0] )*( ( X[2][1] - X[0][1] )*( X[3][2] - X[0][2] ) - ( X[3][1] - X[0][1] )*( X[2][2] - X[0][2] ) )
576  + ( X[1][1] - X[0][1] )*( ( X[3][0] - X[0][0] )*( X[2][2] - X[0][2] ) - ( X[2][0] - X[0][0] )*( X[3][2] - X[0][2] ) )
577  + ( X[1][2] - X[0][2] )*( ( X[2][0] - X[0][0] )*( X[3][1] - X[0][1] ) - ( X[3][0] - X[0][0] )*( X[2][1] - X[0][1] ) );
578 }
579 
580 //*************************************************************************************************
581 
582 template< typename NUM_Q_POINTS >
584 inline
585 void
587 calcN( real64 const (&coords)[3],
588  real64 (& N)[numNodes] )
589 {
590  real64 const r = coords[0];
591  real64 const s = coords[1];
592  real64 const t = coords[2];
593 
594  // single quadrature point (centroid), i.e. r = s = t = 1/4
595  N[0] = 1 - r - s - t;
596  N[1] = r;
597  N[2] = s;
598  N[3] = t;
599 }
600 
601 
602 template< typename NUM_Q_POINTS >
604 inline
605 void
607 calcN( localIndex const q,
608  real64 (& N)[numNodes] )
609 {
610  real64 const pointCoord[3] = {quadratureParentCoords0( q ),
611  quadratureParentCoords1( q ),
612  quadratureParentCoords2( q )};
613 
614  calcN( pointCoord, N );
615 }
616 
617 template< typename NUM_Q_POINTS >
619 inline
621 calcN( localIndex const q,
622  StackVariables const & GEOS_UNUSED_PARAM( stack ),
623  real64 ( & N )[numNodes] )
624 {
625  return calcN( q, N );
626 }
627 
628 //*************************************************************************************************
629 template< typename NUM_Q_POINTS >
632 real64
634  real64 const (&X)[numNodes][3],
635  real64 (& J)[3][3] )
636 {
637  for( int i = 0; i < 3; ++i )
638  {
639  for( int j = 0; j < 3; ++j )
640  {
641  J[i][j] = 0;
642  }
643  }
644  J[0][0] = X[1][0] - X[0][0];
645  J[0][1] = X[2][0] - X[0][0];
646  J[0][2] = X[3][0] - X[0][0];
647 
648  J[1][0] = X[1][1] - X[0][1];
649  J[1][1] = X[2][1] - X[0][1];
650  J[1][2] = X[3][1] - X[0][1];
651 
652  J[2][0] = X[1][2] - X[0][2];
653  J[2][1] = X[2][2] - X[0][2];
654  J[2][2] = X[3][2] - X[0][2];
655  real64 const detJ = LvArray::tensorOps::determinant< 3 >( J );
656  return detJ * quadratureWeight( q );
657 }
658 
659 template< typename NUM_Q_POINTS >
662 real64
664  real64 const (&X)[numNodes][3],
665  StackVariables const & GEOS_UNUSED_PARAM( stack ),
666  real64 (& J)[3][3] )
667 {
668  return calcJacobian( q, X, J );
669 }
670 
671 
672 
673 template< typename NUM_Q_POINTS >
675 inline
676 real64
678 calcGradN( localIndex const q,
679  real64 const (&X)[numNodes][3],
680  real64 (& gradN)[numNodes][3] )
681 {
682 
683  gradN[0][0] = X[1][1]*( X[3][2] - X[2][2] ) - X[2][1]*( X[3][2] - X[1][2] ) + X[3][1]*( X[2][2] - X[1][2] );
684  gradN[0][1] = -X[1][0]*( X[3][2] - X[2][2] ) + X[2][0]*( X[3][2] - X[1][2] ) - X[3][0]*( X[2][2] - X[1][2] );
685  gradN[0][2] = X[1][0]*( X[3][1] - X[2][1] ) - X[2][0]*( X[3][1] - X[1][1] ) + X[3][0]*( X[2][1] - X[1][1] );
686 
687  gradN[1][0] = -X[0][1]*( X[3][2] - X[2][2] ) + X[2][1]*( X[3][2] - X[0][2] ) - X[3][1]*( X[2][2] - X[0][2] );
688  gradN[1][1] = X[0][0]*( X[3][2] - X[2][2] ) - X[2][0]*( X[3][2] - X[0][2] ) + X[3][0]*( X[2][2] - X[0][2] );
689  gradN[1][2] = -X[0][0]*( X[3][1] - X[2][1] ) + X[2][0]*( X[3][1] - X[0][1] ) - X[3][0]*( X[2][1] - X[0][1] );
690 
691  gradN[2][0] = X[0][1]*( X[3][2] - X[1][2] ) - X[1][1]*( X[3][2] - X[0][2] ) + X[3][1]*( X[1][2] - X[0][2] );
692  gradN[2][1] = -X[0][0]*( X[3][2] - X[1][2] ) + X[1][0]*( X[3][2] - X[0][2] ) - X[3][0]*( X[1][2] - X[0][2] );
693  gradN[2][2] = X[0][0]*( X[3][1] - X[1][1] ) - X[1][0]*( X[3][1] - X[0][1] ) + X[3][0]*( X[1][1] - X[0][1] );
694 
695  gradN[3][0] = -X[0][1]*( X[2][2] - X[1][2] ) + X[1][1]*( X[2][2] - X[0][2] ) - X[2][1]*( X[1][2] - X[0][2] );
696  gradN[3][1] = X[0][0]*( X[2][2] - X[1][2] ) - X[1][0]*( X[2][2] - X[0][2] ) + X[2][0]*( X[1][2] - X[0][2] );
697  gradN[3][2] = -X[0][0]*( X[2][1] - X[1][1] ) + X[1][0]*( X[2][1] - X[0][1] ) - X[2][0]*( X[1][1] - X[0][1] );
698 
699  real64 detJ = determinantJacobianTransformation( X );
700  real64 factor = 1.0 / ( detJ );
701 
702  for( int i = 0; i < numNodes; ++i )
703  {
704  for( int j = 0; j < 3; ++j )
705  {
706  gradN[i][j] *= factor;
707  }
708  }
709 
710  return detJ * weight * quadratureWeight( q );
711 }
712 
713 template< typename NUM_Q_POINTS >
715 inline
717 calcGradN( localIndex const q,
718  real64 const (&X)[numNodes][3],
719  StackVariables const & GEOS_UNUSED_PARAM( stack ),
720  real64 ( & gradN )[numNodes][3] )
721 {
722  return calcGradN( q, X, gradN );
723 }
724 
725 template< typename NUM_Q_POINTS >
727 inline
728 real64
730  real64 const (&X)[numNodes][3],
731  real64 (& gradN)[numFaces][3] )
732 {
733 
734  //real64 detJ = determinantJacobianTransformation( X );
735 
736  real64 J[3][3] = {{0}};
737 
738  J[0][0] = X[1][0] - X[0][0];
739  J[0][1] = X[2][0] - X[0][0];
740  J[0][2] = X[3][0] - X[0][0];
741 
742  J[1][0] = X[1][1] - X[0][1];
743  J[1][1] = X[2][1] - X[0][1];
744  J[1][2] = X[3][1] - X[0][1];
745 
746  J[2][0] = X[1][2] - X[0][2];
747  J[2][1] = X[2][2] - X[0][2];
748  J[2][2] = X[3][2] - X[0][2];
749 
750  real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
751 
752  real64 dNdXi[numFaces][3] = {{0}};
753 
754  real64 const r = quadratureParentCoords0( q );
755  real64 const s = quadratureParentCoords1( q );
756  real64 const t = quadratureParentCoords2( q );
757 
758  dNdXi[0][0] = ( 1 - 2 * r - s - t ) * t; // dN0/dr
759  dNdXi[0][1] = -r * t; // dN0/ds
760  dNdXi[0][2] = ( 1 - r - s - 2 * t ) * r; // dN0/dt
761 
762  dNdXi[1][0] = ( 1 - 2 * r - s - t ) * s; // dN1/dr
763  dNdXi[1][1] = ( 1 - r - 2 * s - t ) * r; // dN1/ds
764  dNdXi[1][2] = -r * s; // dN1/dt
765 
766  dNdXi[2][0] = -t * s; // dN2/dr
767  dNdXi[2][1] = ( 1 - r - 2 * s - t ) * t; // dN2/ds
768  dNdXi[2][2] = ( 1 - r - s - 2 * t ) * s; // dN2/dt
769 
770  dNdXi[3][0] = t * s; // dN3/dr
771  dNdXi[3][1] = r * t; // dN3/ds
772  dNdXi[3][2] = r * s; // dN3/dt
773 
774  for( int fi=0; fi<numFaces; ++fi )
775  {
776  gradN[fi][0] = dNdXi[fi][0] * J[0][0] + dNdXi[fi][1] * J[1][0] + dNdXi[fi][2] * J[2][0];
777  gradN[fi][1] = dNdXi[fi][0] * J[0][1] + dNdXi[fi][1] * J[1][1] + dNdXi[fi][2] * J[2][1];
778  gradN[fi][2] = dNdXi[fi][0] * J[0][2] + dNdXi[fi][1] * J[1][2] + dNdXi[fi][2] * J[2][2];
779  }
780 
781  return detJ * weight * quadratureWeight( q );
782 }
783 
784 //*************************************************************************************************
785 
786 template< typename NUM_Q_POINTS >
788 inline
789 real64
792  real64 const (&X)[numNodes][3] )
793 {
794 
795  real64 detJ = determinantJacobianTransformation( X );
796 
797  return detJ * weight * quadratureWeight( q );
798 }
799 
801 
802 
804 template< typename NUM_Q_POINTS >
806  public FiniteElementBase
807 {
808 public:
809 
812 
815 
821  {}
822 
824  virtual ~H1_Tetrahedron_Lagrange1_Gauss() override final = default;
825 
831  {
832  return static_cast< H1_Tetrahedron_Lagrange1_Gauss_impl< NUM_Q_POINTS > * >(this);
833  }
834 
840  {
841  return static_cast< H1_Tetrahedron_Lagrange1_Gauss_impl< NUM_Q_POINTS > const * >(this);
842  }
843 
844 };
845 
852 
859 }
860 }
861 
862 #endif //GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1TETRAHEDRONLAGRANGE1GAUSS_HPP_
#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_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 localIndex numSupportPoints
The number of support points per element.
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 real64 calcGradFaceBubbleN(localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numFaces][3])
Calculate the bubble function 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 calcFaceBubbleN(localIndex const q, real64(&N)[numFaces])
Calculate face bubble functions values for each face at a quadrature point.
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.
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 localIndex getMaxSupportPoints()
Get the maximum number of support 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.
constexpr static int numSamplingPointsPerDirection
Number of sampling points per direction.
constexpr static localIndex numFaces
Number of faces in the 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.
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 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3])
Calculate the integration weights for 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.
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.
constexpr static localIndex numSupportPoints
Number of support points in the element.
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], StackVariables const &stack, real64(&gradN)[numNodes][3])
Calculate the shape functions derivatives wrt the physical coordinates.
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.
constexpr static int numSamplingPoints
Total number of sampling points.
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.
typename Base::template MeshData< SubregionType > MeshData
Mesh data structure for the element.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
H1_Tetrahedron_Lagrange1_Gauss_impl< NUM_Q_POINTS > * getImpl()
Get the device-compatible implementation type.
H1_Tetrahedron_Lagrange1_Gauss_impl< NUM_Q_POINTS > const * getImpl() const
Get the device-compatible implementation type.
virtual ~H1_Tetrahedron_Lagrange1_Gauss() override final=default
Destructor.
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.