GEOS
H1_QuadrilateralFace_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_H1QUADRILATERALFACELAGRANGE1GAUSSLEGENDRE2_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1QUADRILATERALFACELAGRANGE1GAUSSLEGENDRE2_HPP_
22 
23 #include "FiniteElementBase.hpp"
24 #include "LagrangeBasis1.hpp"
25 
26 
27 namespace geos
28 {
29 namespace finiteElement
30 {
31 
51 {
52 public:
53 
56 
58  constexpr static localIndex numNodes = 4;
60  constexpr static localIndex maxSupportPoints = numNodes;
61 
63  constexpr static localIndex numQuadraturePoints = 4;
64 
67  {}
68 
70  virtual localIndex getNumQuadraturePoints() const override
71  {
72  return numQuadraturePoints;
73  }
74 
82  {
83  GEOS_UNUSED_VAR( stack );
84  return numQuadraturePoints;
85  }
86 
88  virtual localIndex getNumSupportPoints() const override
89  {
90  return numNodes;
91  }
92 
94  virtual localIndex getMaxSupportPoints() const override
95  {
96  return maxSupportPoints;
97  }
98 
106  {
107  GEOS_UNUSED_VAR( stack );
108  return numNodes;
109  }
110 
117  static void calcN( real64 const (&coords)[2],
118  real64 ( &N )[numNodes] );
119 
128  static void calcN( localIndex const q,
129  real64 ( &N )[numNodes] );
130 
140  inline
141  static void calcN( localIndex const q,
142  StackVariables const & stack,
143  real64 ( &N )[numNodes] );
144 
152  static void calcBubbleN( real64 const (&pointCoord)[2],
153  real64 (& N)[1] )
154  {
156  }
157 
165  inline
166  static void calcBubbleN( localIndex const q,
167  real64 (& N)[1] )
168  {
169  real64 const qCoords[2] = { quadratureFactor *parentCoords0( q ),
170  quadratureFactor *parentCoords1( q ) };
171 
172  calcBubbleN( qCoords, N );
173  }
174 
184  real64 const (&X)[numNodes][3] );
185 
195  template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, bool UPPER >
197  inline
198  static void addGradGradStabilization( StackVariables const & stack,
199  real64 ( &matrix )
200  [maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT]
201  [maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT],
202  real64 const & scaleFactor );
203 
210  inline
211  static void getPermutation( int (& permutation)[numNodes] )
212  {
213  permutation[0] = 0;
214  permutation[1] = 1;
215  permutation[2] = 3;
216  permutation[3] = 2;
217  }
218 
219 private:
221  constexpr static real64 parentArea = 4.0;
222 
224  constexpr static real64 weight = parentArea / numQuadraturePoints;
225 
229  constexpr static real64 quadratureFactor = 1.0 / 1.732050807568877293528;
230 
238  template< typename T >
240  inline
241  constexpr static T linearMap( T const i, T const j )
242  {
243  return i + 2 * j;
244  }
245 
253  inline
254  constexpr static real64 parentCoords0( localIndex const a )
255  {
256  return -1.0 + 2.0 * (a & 1);
257  }
258 
266  inline
267  constexpr static real64 parentCoords1( localIndex const a )
268  {
269  return -1.0 + ( a & 2 );
270  }
271 
272 };
273 
275 
276 template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, bool UPPER >
278 inline
280  addGradGradStabilization( StackVariables const & stack,
281  real64 ( & matrix )
282  [maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT]
283  [maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT],
284  real64 const & scaleFactor )
285 {
286  GEOS_UNUSED_VAR( stack );
287  GEOS_UNUSED_VAR( matrix );
288  GEOS_UNUSED_VAR( scaleFactor );
289 }
290 
292 inline
293 void
295  calcN( real64 const (&coords)[2],
296  real64 (& N)[numNodes] )
297 {
298  for( localIndex a=0; a<numNodes; ++a )
299  {
300  N[a] = 0.25 *
301  ( 1 + quadratureFactor*coords[0]*parentCoords0( a ) ) *
302  ( 1 + quadratureFactor*coords[1]*parentCoords1( a ) );
303  }
304 }
306 inline
307 void
309  calcN( localIndex const q,
310  real64 (& N)[numNodes] )
311 {
312  for( localIndex a=0; a<numNodes; ++a )
313  {
314  N[a] = 0.25 *
315  ( 1 + quadratureFactor*parentCoords0( q )*parentCoords0( a ) ) *
316  ( 1 + quadratureFactor*parentCoords1( q )*parentCoords1( a ) );
317  }
318 }
319 
321 inline
323  calcN( localIndex const q,
324  StackVariables const & GEOS_UNUSED_PARAM( stack ),
325  real64 ( & N )[numNodes] )
326 {
327  return calcN( q, N );
328 }
329 
330 //*************************************************************************************************
331 
333 inline
334 real64
337  real64 const (&X)[numNodes][3] )
338 {
339  real64 dXdXi[3][2] = {{0}};
340 
341  real64 const quadratureCoords[2] = { quadratureFactor *parentCoords0( q ),
342  quadratureFactor *parentCoords1( q ) };
343 
344  real64 const psi0[2] = { 0.5*( 1.0 - quadratureCoords[0] ),
345  0.5*( 1.0 + quadratureCoords[0] ) };
346  real64 const psi1[2] = { 0.5*( 1.0 - quadratureCoords[1] ),
347  0.5*( 1.0 + quadratureCoords[1] ) };
348  constexpr real64 dpsi[2] = { -0.5, 0.5 };
349 
350  for( int a=0; a<2; ++a )
351  {
352  for( int b=0; b<2; ++b )
353  {
354  real64 const dNdXi[2] = { dpsi[a] * psi1[b],
355  psi0[a] * dpsi[b] };
356 
357  localIndex const nodeIndex = linearMap( a, b );
358 
359  for( int i = 0; i < 3; ++i )
360  {
361  for( int j = 0; j < 2; ++j )
362  {
363  dXdXi[i][j] = dXdXi[i][j] + dNdXi[ j ] * X[nodeIndex][i];
364  }
365  }
366  }
367  }
368 
369  real64 const detJ = pow( dXdXi[1][0] * dXdXi[2][1] - dXdXi[2][0] * dXdXi[1][1], 2.0 )
370  + pow( dXdXi[2][0] * dXdXi[0][1] - dXdXi[0][0] * dXdXi[2][1], 2.0 )
371  + pow( dXdXi[0][0] * dXdXi[1][1] - dXdXi[1][0] * dXdXi[0][1], 2.0 );
372 
373  return sqrt( detJ ) * weight;
374 }
375 
377 
378 }
379 }
380 #endif //GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1QUADRILATERALFACELAGRANGE1GAUSSLEGENDRE2_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
Base class for FEM element implementations.
constexpr static localIndex maxSupportPoints
The maximum number of support points per element.
static GEOS_HOST_DEVICE void getPermutation(int(&permutation)[numNodes])
Calculate the node permutation between the parent element and the geometric element....
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
constexpr static localIndex numNodes
The number of nodes/support 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.
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 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 calcBubbleN(real64 const (&pointCoord)[2], real64(&N)[1])
Calculate shape bubble functions values at a given point in the parent space.
virtual GEOS_HOST_DEVICE localIndex getMaxSupportPoints() const override
Get the maximum number of support points for this element.
static GEOS_HOST_DEVICE void calcBubbleN(localIndex const q, real64(&N)[1])
Calculate shape bubble functions values at a quadrature point.
constexpr static localIndex numQuadraturePoints
The number of quadrature points per element.
static GEOS_HOST_DEVICE void calcN(real64 const (&coords)[2], real64(&N)[numNodes])
Calculate shape functions values at a single point.
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 addGradGradStabilization(StackVariables const &stack, real64(&matrix)[maxSupportPoints *NUMDOFSPERTRIALSUPPORTPOINT][maxSupportPoints *NUMDOFSPERTRIALSUPPORTPOINT], real64 const &scaleFactor)
Empty method, here for compatibility with methods that require a stabilization of the grad-grad bilin...
virtual GEOS_HOST_DEVICE localIndex getNumQuadraturePoints() const override
Virtual getter for the number of quadrature points per element.
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
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void valueBubble(real64 const (&coords)[2], real64(&N)[1])
The value of the bubble basis function evaluated at a point along the axes.