GEOS
H1_Pyramid_Lagrange1_Gauss5.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_H1PYRAMIDLAGRANGE1GAUSS5_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1PYRAMIDLAGRANGE1GAUSS5_HPP_
22 
23 #include "FiniteElementBase.hpp"
24 #include "LagrangeBasis1.hpp"
25 
26 
27 namespace geos
28 {
29 namespace finiteElement
30 {
31 
59 {
60 public:
61 
64 
66  constexpr static localIndex numNodes = 5;
67 
69  constexpr static localIndex numFaces = 5;
70 
72  constexpr static localIndex maxSupportPoints = numNodes;
73 
75  constexpr static localIndex numQuadraturePoints = 5;
76 
79 
81  virtual ~H1_Pyramid_Lagrange1_Gauss5() override
82  {}
83 
85  virtual localIndex getNumQuadraturePoints() const override
86  {
87  return numQuadraturePoints;
88  }
89 
97  {
98  GEOS_UNUSED_VAR( stack );
99  return numQuadraturePoints;
100  }
101 
103  virtual localIndex getNumSupportPoints() const override
104  {
105  return numNodes;
106  }
107 
115  {
116  GEOS_UNUSED_VAR( stack );
117  return numNodes;
118  }
119 
121  virtual localIndex getMaxSupportPoints() const override
122  {
123  return maxSupportPoints;
124  }
125 
134  static void getSamplingPointCoordInParentSpace( int const & linearIndex,
135  real64 (& samplingPointCoord)[3] )
136  {
137  GEOS_UNUSED_VAR( linearIndex, samplingPointCoord );
138  GEOS_ERROR( " Element type not supported." );
139  }
140 
148  static void calcN( real64 const (&pointCoord)[3],
149  real64 ( &N )[numNodes] );
150 
159  static void calcN( localIndex const q,
160  real64 ( &N )[numNodes] );
161 
171  inline
172  static void calcN( localIndex const q,
173  StackVariables const & stack,
174  real64 ( &N )[numNodes] );
175 
185  static void calcFaceBubbleN( real64 const (&pointCoord)[3],
186  real64 (& N)[numFaces] )
187  {
188  GEOS_UNUSED_VAR( pointCoord, N );
189  GEOS_ERROR( "Unsupported bubble functions for pyramid elements" );
190  }
191 
200  inline
201  static void calcFaceBubbleN( localIndex const q,
202  real64 (& N)[numFaces] )
203  {
204  GEOS_UNUSED_VAR( q, N );
205  GEOS_ERROR( "Unsupported bubble functions for pyramid elements" );
206  }
207 
218  static real64 calcGradN( localIndex const q,
219  real64 const (&X)[numNodes][3],
220  real64 ( &gradN )[numNodes][3] );
221 
233  inline
234  static real64 calcGradN( localIndex const q,
235  real64 const (&X)[numNodes][3],
236  StackVariables const & stack,
237  real64 ( &gradN )[numNodes][3] );
238 
250  real64 const (&X)[numNodes][3],
251  real64 ( &gradN )[numFaces][3] );
252 
262  real64 const (&X)[numNodes][3] );
263 
274  real64 const (&X)[numNodes][3],
275  StackVariables const & stack )
276  { GEOS_UNUSED_VAR( stack ); return transformedQuadratureWeight( q, X ); }
277 
287  static real64 invJacobianTransformation( int const q,
288  real64 const (&X)[numNodes][3],
289  real64 ( & J )[3][3] )
290  {
291  jacobianTransformation( q, X, J );
292  return LvArray::tensorOps::invert< 3 >( J );
293  }
294 
295 private:
297  constexpr static real64 weight = 81.0 / 100.0;
298 
300  constexpr static real64 weightDelta = 125.0 / 27.0 - weight;
301 
305  constexpr static real64 quadratureCrossSectionCoord = 0.584237394672177;
306 
309  constexpr static real64 quadratureLongitudinalCoordNeg = -2.0 / 3.0;
310 
313  constexpr static real64 quadratureLongitudinalCoordDelta = 16.0 / 15.0;
314 
322  template< typename T >
324  inline
325  constexpr static T linearMap( T const i, T const j )
326  {
327  return i + 2 * j;
328  }
329 
337  inline
338  constexpr static real64 parentCoords0( localIndex const a )
339  {
340  return -1.0 + 2.0 * ( a & 1 ) + 0.25 * ( a & 4 );
341  }
342 
350  inline
351  constexpr static real64 parentCoords1( localIndex const a )
352  {
353  return -1.0 + ( a & 2 ) + 0.25 * ( a & 4 );
354  }
355 
363  inline
364  constexpr static real64 parentCoords2( localIndex const a )
365  {
366  return -1.0 + 0.5 * ( a & 4 );
367  }
368 
376  inline
377  constexpr static real64 quadratureParentCoords0( localIndex const q )
378  {
379  return parentCoords0( q ) * quadratureCrossSectionCoord;
380  }
381 
389  inline
390  constexpr static real64 quadratureParentCoords1( localIndex const q )
391  {
392  return parentCoords1( q ) * quadratureCrossSectionCoord;
393  }
394 
402  inline
403  constexpr static real64 quadratureParentCoords2( localIndex const q )
404  {
405  return quadratureLongitudinalCoordNeg + 0.5 * ( 1 + parentCoords2( q ) ) * quadratureLongitudinalCoordDelta;
406  }
407 
414  inline
415  constexpr static real64 quadratureWeight( localIndex const q )
416  {
417  return weight + 0.5 * ( 1 + parentCoords2( q ) ) * weightDelta;
418  }
419 
427  static void jacobianTransformation( int const q,
428  real64 const (&X)[numNodes][3],
429  real64 ( &J )[3][3] );
430 
441  static void
442  applyJacobianTransformationToShapeFunctionsDerivatives( int const q,
443  real64 const ( &invJ )[3][3],
444  real64 ( &gradN )[numNodes][3] );
445 
446 };
447 
449 
451 inline
452 void
453 H1_Pyramid_Lagrange1_Gauss5::
454  jacobianTransformation( int const q,
455  real64 const (&X)[numNodes][3],
456  real64 ( & J )[3][3] )
457 {
458  real64 const quadratureCoords[3] = { quadratureParentCoords0( q ),
459  quadratureParentCoords1( q ),
460  quadratureParentCoords2( q ) };
461 
462  real64 const psi0[2] = { 0.5 - 0.5*quadratureCoords[0],
463  0.5 + 0.5*quadratureCoords[0] };
464  real64 const psi1[2] = { 0.5 - 0.5*quadratureCoords[1],
465  0.5 + 0.5*quadratureCoords[1] };
466  real64 const psi2 = 0.5 - 0.5*quadratureCoords[2];
467  constexpr real64 dpsi[2] = { -0.5, 0.5 };
468 
469  // Contributions from basis functions paired with base nodes
470  for( localIndex a=0; a<2; ++a )
471  {
472  for( localIndex b=0; b<2; ++b )
473  {
474  real64 const dNdXi[3] = { dpsi[a] * psi1[b] * psi2,
475  psi0[a] * dpsi[b] * psi2,
476  psi0[a] * psi1[b] * dpsi[0] };
477  localIndex const nodeIndex = linearMap( a, b );
478  for( int i = 0; i < 3; ++i )
479  {
480  for( int j = 0; j < 3; ++j )
481  {
482  J[i][j] = J[i][j] + dNdXi[ j ] * X[nodeIndex][i];
483  }
484  }
485  }
486  }
487 
488  // Contribution from the basis function paired with the apex nodes
489  for( int i = 0; i < 3; ++i )
490  {
491  J[i][2] = J[i][2] + dpsi[1] * X[4][i];
492  }
493 }
494 
495 //*************************************************************************************************
496 
498 inline
499 void
500 H1_Pyramid_Lagrange1_Gauss5::
501  applyJacobianTransformationToShapeFunctionsDerivatives( int const q,
502  real64 const ( &invJ )[3][3],
503  real64 (& gradN)[numNodes][3] )
504 {
505  real64 const quadratureCoords[3] = { quadratureParentCoords0( q ),
506  quadratureParentCoords1( q ),
507  quadratureParentCoords2( q ) };
508 
509  real64 const psi0[2] = { 0.5*( 1.0 - quadratureCoords[0] ),
510  0.5*( 1.0 + quadratureCoords[0] ) };
511  real64 const psi1[2] = { 0.5*( 1.0 - quadratureCoords[1] ),
512  0.5*( 1.0 + quadratureCoords[1] ) };
513  real64 const psi2 = 0.5*( 1.0 - quadratureCoords[2]);
514  constexpr real64 dpsi[2] = { -0.5, 0.5 };
515 
516  // Contributions from basis functions paired with base nodes
517  for( localIndex a=0; a<2; ++a )
518  {
519  for( localIndex b=0; b<2; ++b )
520  {
521  real64 const dNdXi[3] = { dpsi[a] * psi1[b] * psi2,
522  psi0[a] * dpsi[b] * psi2,
523  psi0[a] * psi1[b] * dpsi[0] };
524  localIndex const nodeIndex = linearMap( a, b );
525  for( int i = 0; i < 3; ++i )
526  {
527  gradN[nodeIndex][i] = 0.0;
528  for( int j = 0; j < 3; ++j )
529  {
530  gradN[nodeIndex][i] = gradN[nodeIndex][i] + dNdXi[ j ] * invJ[j][i];
531  }
532  }
533  }
534  }
535 
536  // Contribution from the basis function paired with the apex nodes
537  for( int i = 0; i < 3; ++i )
538  {
539  gradN[4][i] = dpsi[1] * invJ[2][i];
540  }
541 }
542 
543 //*************************************************************************************************
544 
546 inline
547 void
549  calcN( real64 const ( &pointCoord )[3],
550  real64 ( & N )[numNodes] )
551 {
552  N[0] = 0.125*( 1.0 - pointCoord[0] ) * ( 1.0 - pointCoord[1] ) * ( 1.0 - pointCoord[2] );
553  N[1] = 0.125*( 1.0 + pointCoord[0] ) * ( 1.0 - pointCoord[1] ) * ( 1.0 - pointCoord[2] );
554  N[2] = 0.125*( 1.0 - pointCoord[0] ) * ( 1.0 + pointCoord[1] ) * ( 1.0 - pointCoord[2] );
555  N[3] = 0.125*( 1.0 + pointCoord[0] ) * ( 1.0 + pointCoord[1] ) * ( 1.0 - pointCoord[2] );
556  N[4] = 0.5*( 1.0 + pointCoord[2] );
557 }
558 
560 inline
561 void
563  calcN( localIndex const q,
564  real64 ( & N )[numNodes] )
565 {
566  real64 const xi[3] = { quadratureParentCoords0( q ),
567  quadratureParentCoords1( q ),
568  quadratureParentCoords2( q ) };
569 
570  calcN( xi, N );
571 }
572 
574 inline
576  calcN( localIndex const q,
577  StackVariables const & GEOS_UNUSED_PARAM( stack ),
578  real64 ( & N )[numNodes] )
579 {
580  return calcN( q, N );
581 }
582 
583 //*************************************************************************************************
584 
586 inline
588  real64 const (&X)[numNodes][3],
589  real64 (& gradN)[numNodes][3] )
590 {
591  real64 J[3][3] = {{0}};
592 
593  jacobianTransformation( q, X, J );
594 
595  real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
596 
597  applyJacobianTransformationToShapeFunctionsDerivatives( q, J, gradN );
598 
599  return detJ * quadratureWeight( q );
600 }
601 
603 inline
605  calcGradN( localIndex const q,
606  real64 const (&X)[numNodes][3],
607  StackVariables const & GEOS_UNUSED_PARAM( stack ),
608  real64 ( & gradN )[numNodes][3] )
609 {
610  return calcGradN( q, X, gradN );
611 }
612 
614 inline
615 real64
617  real64 const (&X)[numNodes][3],
618  real64 (& gradN)[numFaces][3] )
619 {
620  GEOS_UNUSED_VAR( q, X, gradN );
621  GEOS_ERROR( "Unsupported bubble functions for pyramid elements" );
622  return 0.0;
623 }
624 
625 //*************************************************************************************************
626 
628 inline
629 real64
632  real64 const (&X)[numNodes][3] )
633 {
634  real64 J[3][3] = {{0}};
635 
636  jacobianTransformation( q, X, J );
637 
638  return LvArray::tensorOps::determinant< 3 >( J ) * quadratureWeight( q );
639 }
640 
642 
643 }
644 }
645 #endif //GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1PYRAMIDLAGRANGE1GAUSS5_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
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
Definition: Logger.hpp:157
Base class for FEM element implementations.
constexpr static int numSamplingPointsPerDirection
Number of sampling points.
virtual GEOS_HOST_DEVICE localIndex getNumSupportPoints() const override
Virtual getter for the number of support points per element.
constexpr static localIndex numFaces
The number of faces/support points per element.
static GEOS_HOST_DEVICE real64 calcGradFaceBubbleN(localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numFaces][3])
Calculate the shape bubble function derivatives wrt the physical coordinates.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
virtual GEOS_HOST_DEVICE localIndex getNumQuadraturePoints() const override
Virtual getter for the number of quadrature 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.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints(StackVariables const &stack)
Get the number of support points.
constexpr static int numSamplingPoints
The number of sampling points per 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 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 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 calcN(real64 const (&pointCoord)[3], real64(&N)[numNodes])
Calculate shape functions values at a single point.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcFaceBubbleN(real64 const (&pointCoord)[3], real64(&N)[numFaces])
Calculate shape functions values for each support face at a given point in the parent space.
constexpr static localIndex numNodes
The number of nodes/support points per element.
virtual GEOS_HOST_DEVICE localIndex getMaxSupportPoints() const override
Get the maximum number of support points for this element.
constexpr static localIndex maxSupportPoints
The maximum number of support points per element.
constexpr static localIndex numQuadraturePoints
The number of quadrature points per element.
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...
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.
static GEOS_HOST_DEVICE void calcFaceBubbleN(localIndex const q, real64(&N)[numFaces])
Calculate shape functions values for each support face at a quadrature point.
static GEOS_HOST_DEVICE real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3])
Calculate the integration weights for a quadrature 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