GEOS
H1_Wedge_Lagrange1_Gauss6.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_H1WEDGELAGRANGE1GAUSS6_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1WEDGELAGRANGE1GAUSS6_HPP_
22 
23 #include "FiniteElementBase.hpp"
24 #include "LagrangeBasis1.hpp"
25 
26 
27 namespace geos
28 {
29 namespace finiteElement
30 {
31 
56 {
57 public:
58 
61 
63  constexpr static localIndex numNodes = 6;
64 
66  constexpr static localIndex numFaces = 5;
67 
69  constexpr static localIndex maxSupportPoints = numNodes;
70 
72  constexpr static localIndex numQuadraturePoints = 6;
73 
76 
78  virtual ~H1_Wedge_Lagrange1_Gauss6() override
79  {}
80 
82  virtual localIndex getNumQuadraturePoints() const override
83  {
84  return numQuadraturePoints;
85  }
86 
94  {
95  GEOS_UNUSED_VAR( stack );
96  return numQuadraturePoints;
97  }
98 
100  virtual localIndex getNumSupportPoints() const override
101  {
102  return numNodes;
103  }
104 
106  virtual localIndex getMaxSupportPoints() const override
107  {
108  return maxSupportPoints;
109  }
110 
118  {
119  GEOS_UNUSED_VAR( stack );
120  return numNodes;
121  }
122 
131  static void getSamplingPointCoordInParentSpace( int const & linearIndex,
132  real64 (& samplingPointCoord)[3] )
133  {
134  int const i0 = linearIndex % numSamplingPointsPerDirection;
135  int const i1 = ( (linearIndex - i0)/numSamplingPointsPerDirection ) % numSamplingPointsPerDirection;
136  int const i2 = ( (linearIndex - i0)/numSamplingPointsPerDirection - i1 ) / numSamplingPointsPerDirection;
137 
138  real64 const step = 1. / ( numSamplingPointsPerDirection - 1 );
139 
140  real64 const r = i0 * step;
141  real64 const s = i1 * step;
142  real64 const t = i2 * 2 * step;
143  if( (r+s) <= 1 )
144  {
145  samplingPointCoord[0] = r;
146  samplingPointCoord[1] = s;
147  samplingPointCoord[2] = t;
148  }
149  else
150  {
151  // if outside of the triangle need to reproject it. Points will be doubled though.
152  samplingPointCoord[0] = 1 - r;
153  samplingPointCoord[1] = 1 - s;
154  samplingPointCoord[2] = t;
155  }
156  }
157 
164  inline
165  static void calcN( real64 const (&coords)[3],
166  real64 ( &N )[numNodes] );
167 
168 
177  static void calcN( localIndex const q,
178  real64 ( &N )[numNodes] );
179 
189  inline
190  static void calcN( localIndex const q,
191  StackVariables const & stack,
192  real64 ( &N )[numNodes] );
193 
203  static void calcFaceBubbleN( real64 const (&pointCoord)[3],
204  real64 (& N)[numFaces] )
205  {
206  GEOS_UNUSED_VAR( pointCoord, N );
207  GEOS_ERROR( "Unsupported bubble functions for wedge elements" );
208  }
209 
218  inline
219  static void calcFaceBubbleN( localIndex const q,
220  real64 (& N)[numFaces] )
221  {
222  GEOS_UNUSED_VAR( q, N );
223  GEOS_ERROR( "Unsupported bubble functions for wedge elements" );
224  }
225 
236  static real64 calcGradN( localIndex const q,
237  real64 const (&X)[numNodes][3],
238  real64 ( &gradN )[numNodes][3] );
239 
251  inline
252  static real64 calcGradN( localIndex const q,
253  real64 const (&X)[numNodes][3],
254  StackVariables const & stack,
255  real64 ( &gradN )[numNodes][3] );
256 
268  real64 const (&X)[numNodes][3],
269  real64 ( &gradN )[numFaces][3] );
270 
280  real64 const (&X)[numNodes][3] );
281 
292  real64 const (&X)[numNodes][3],
293  StackVariables const & stack )
294  { GEOS_UNUSED_VAR( stack ); return transformedQuadratureWeight( q, X ); }
295 
305  static real64 invJacobianTransformation( int const q,
306  real64 const (&X)[numNodes][3],
307  real64 ( & J )[3][3] )
308  {
309  jacobianTransformation( q, X, J );
310  return LvArray::tensorOps::invert< 3 >( J );
311  }
312 
313 
314 private:
316  constexpr static real64 parentVolume = 1.0;
317 
319  constexpr static real64 weight = parentVolume / numQuadraturePoints;
320 
323  constexpr static real64 quadratureCrossSectionCoord = 1.0 / 6.0;
324 
328  constexpr static real64 quadratureLongitudinalCoord = 1.0 / 1.732050807568877293528;
329 
337  template< typename T >
339  inline
340  constexpr static T linearMap( T const indexT, T const indexL )
341  {
342  return 2 * indexT + indexL;
343  }
344 
352  inline
353  constexpr static real64 parentCoords0( localIndex const a )
354  {
355  return 0.5 * ( a & 2 );
356  }
357 
365  inline
366  constexpr static real64 parentCoords1( localIndex const a )
367  {
368  return 0.25 * ( a & 4 );
369  }
370 
378  inline
379  constexpr static real64 parentCoords2( localIndex const a )
380  {
381  return -1.0 + 2 * ( a & 1 );
382  }
383 
391  inline
392  constexpr static real64 quadratureParentCoords0( localIndex const q )
393  {
394  return quadratureCrossSectionCoord + 0.5 * parentCoords0( q );
395  }
396 
404  inline
405  constexpr static real64 quadratureParentCoords1( localIndex const q )
406  {
407  return quadratureCrossSectionCoord + 0.5 * parentCoords1( q );
408  }
409 
417  inline
418  constexpr static real64 quadratureParentCoords2( localIndex const q )
419  {
420  return parentCoords2( q ) * quadratureLongitudinalCoord;
421  }
422 
430  static void jacobianTransformation( int const q,
431  real64 const (&X)[numNodes][3],
432  real64 ( &J )[3][3] );
433 
444  static void
445  applyJacobianTransformationToShapeFunctionsDerivatives( int const q,
446  real64 const ( &invJ )[3][3],
447  real64 ( &gradN )[numNodes][3] );
448 
449 };
450 
452 
454 inline
455 void
456 H1_Wedge_Lagrange1_Gauss6::
457  jacobianTransformation( int const q,
458  real64 const (&X)[numNodes][3],
459  real64 ( & J )[3][3] )
460 {
461  real64 const r = quadratureParentCoords0( q );
462  real64 const s = quadratureParentCoords1( q );
463  real64 const xi = quadratureParentCoords2( q );
464 
465  real64 const psiTRI[3] = { 1.0 - r - s, r, s };
466  real64 const psiLIN[2] = { 0.5 - 0.5*xi, 0.5 + 0.5*xi };
467  constexpr real64 dpsiTRI[2][3] = { { -1.0, 1.0, 0.0 }, { -1.0, 0.0, 1.0 } };
468  constexpr real64 dpsiLIN[2] = { -0.5, 0.5 };
469 
470  for( localIndex a=0; a<3; ++a )
471  {
472  for( localIndex b=0; b<2; ++b )
473  {
474  real64 const dNdXi[3] = { dpsiTRI[0][a] * psiLIN[b],
475  dpsiTRI[1][a] * psiLIN[b],
476  psiTRI[a] * dpsiLIN[b] };
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 
489 //*************************************************************************************************
490 
492 inline
493 void
494 H1_Wedge_Lagrange1_Gauss6::
495  applyJacobianTransformationToShapeFunctionsDerivatives( int const q,
496  real64 const ( &invJ )[3][3],
497  real64 (& gradN)[numNodes][3] )
498 {
499  real64 const r = quadratureParentCoords0( q );
500  real64 const s = quadratureParentCoords1( q );
501  real64 const xi = quadratureParentCoords2( q );
502 
503  real64 const psiTRI[3] = { 1.0 - r - s, r, s };
504  real64 const psiLIN[2] = { 0.5 - 0.5*xi, 0.5 + 0.5*xi };
505  constexpr real64 dpsiTRI[2][3] = { { -1.0, 1.0, 0.0 }, { -1.0, 0.0, 1.0 } };
506  constexpr real64 dpsiLIN[2] = { -0.5, 0.5 };
507 
508  for( localIndex a=0; a<3; ++a )
509  {
510  for( localIndex b=0; b<2; ++b )
511  {
512  real64 const dNdXi[3] = { dpsiTRI[0][a] * psiLIN[b],
513  dpsiTRI[1][a] * psiLIN[b],
514  psiTRI[a] * dpsiLIN[b] };
515  localIndex const nodeIndex = linearMap( a, b );
516  for( int i = 0; i < 3; ++i )
517  {
518  gradN[nodeIndex][i] = 0.0;
519  for( int j = 0; j < 3; ++j )
520  {
521  gradN[nodeIndex][i] = gradN[nodeIndex][i] + dNdXi[ j ] * invJ[j][i];
522  }
523  }
524  }
525  }
526 }
527 
528 //*************************************************************************************************
529 
531 inline
532 void
533 H1_Wedge_Lagrange1_Gauss6::calcN( real64 const (&coords)[3],
534  real64 (& N)[numNodes] )
535 {
536  real64 const r = coords[0];
537  real64 const s = coords[1];
538  real64 const xi = coords[2];
539 
540  N[0] = 0.5*( 1.0 - r - s ) * ( 1.0 - xi );
541  N[1] = 0.5*( 1.0 - r - s ) * ( 1.0 + xi );
542  N[2] = 0.5* r * ( 1.0 - xi );
543  N[3] = 0.5* r * ( 1.0 + xi );
544  N[4] = 0.5* s * ( 1.0 - xi );
545  N[5] = 0.5* s * ( 1.0 + xi );
546 
547 }
548 
549 
551 inline
552 void
554  calcN( localIndex const q,
555  real64 (& N)[numNodes] )
556 {
557  real64 const pointCoord[3] = {quadratureParentCoords0( q ),
558  quadratureParentCoords1( q ),
559  quadratureParentCoords2( q )};
560 
561  calcN( pointCoord, N );
562 }
563 
565 inline
567  calcN( localIndex const q,
568  StackVariables const & GEOS_UNUSED_PARAM( stack ),
569  real64 ( & N )[numNodes] )
570 {
571  return calcN( q, N );
572 }
573 
574 //*************************************************************************************************
575 
577 inline
578 real64
580  calcGradN( localIndex const q,
581  real64 const (&X)[numNodes][3],
582  real64 (& gradN)[numNodes][3] )
583 {
584  real64 J[3][3] = {{0}};
585 
586  jacobianTransformation( q, X, J );
587 
588  real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
589 
590  applyJacobianTransformationToShapeFunctionsDerivatives( q, J, gradN );
591 
592  return detJ * weight;
593 }
594 
596 inline
598  calcGradN( localIndex const q,
599  real64 const (&X)[numNodes][3],
600  StackVariables const & GEOS_UNUSED_PARAM( stack ),
601  real64 ( & gradN )[numNodes][3] )
602 {
603  return calcGradN( q, X, gradN );
604 }
605 
607 inline
608 real64
610  real64 const (&X)[numNodes][3],
611  real64 (& gradN)[numFaces][3] )
612 {
613  GEOS_UNUSED_VAR( q, X, gradN );
614  GEOS_ERROR( "Unsupported bubble functions for wedge elements" );
615  return 0.0;
616 }
617 
618 //*************************************************************************************************
619 
621 inline
622 real64
625  real64 const (&X)[numNodes][3] )
626 {
627  real64 J[3][3] = {{0}};
628 
629  jacobianTransformation( q, X, J );
630 
631  return LvArray::tensorOps::determinant< 3 >( J ) * weight;
632 }
633 
635 
636 }
637 }
638 
639 #endif //GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1WEDGELAGRANGE1GAUSS6_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.
static GEOS_HOST_DEVICE void calcFaceBubbleN(localIndex const q, real64(&N)[numFaces])
Calculate shape functions values for each support face at a quadrature point.
virtual GEOS_HOST_DEVICE localIndex getNumSupportPoints() const override
Virtual getter for the number of support points per element.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void getSamplingPointCoordInParentSpace(int const &linearIndex, real64(&samplingPointCoord)[3])
Get the Sampling Point Coord In the Parent Space.
constexpr static localIndex numNodes
The number of nodes/support points per element.
constexpr static localIndex numFaces
The number of faces/support points per element.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
constexpr static localIndex numQuadraturePoints
The number of quadrature 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.
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, StackVariables const &stack, 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 real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3])
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 void calcN(real64 const (&coords)[3], real64(&N)[numNodes])
Calculate shape functions values at a single 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.
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.
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 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 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.
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