GEOS
H1_TriangleFace_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_H1TRIANGLEFACELAGRANGE1GAUSS_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1TRIANGLEFACELAGRANGE1GAUSS_HPP_
22 
23 #include "FiniteElementBase.hpp"
24 #include "LagrangeBasis1.hpp"
25 
26 
27 namespace geos
28 {
29 namespace finiteElement
30 {
31 
50 template< typename NUM_Q_POINTS >
52 {
53 public:
54 
56  static_assert( ( NUM_Q_POINTS::value == 1 ||
57  NUM_Q_POINTS::value == 4 ||
58  NUM_Q_POINTS::value == 6 ),
59  "NUM_Q_POINTS::value must be 1, 4, or 6!" );
60 
63 
65  constexpr static localIndex numNodes = 3;
67  constexpr static localIndex maxSupportPoints = numNodes;
68 
70  constexpr static localIndex numQuadraturePoints = NUM_Q_POINTS::value;
71 
73  virtual ~H1_TriangleFace_Lagrange1_Gauss() override
74  {}
75 
77  virtual localIndex getNumQuadraturePoints() const override
78  {
79  return numQuadraturePoints;
80  }
81 
89  {
90  GEOS_UNUSED_VAR( stack );
91  return numQuadraturePoints;
92  }
93 
95  virtual localIndex getNumSupportPoints() const override
96  {
97  return numNodes;
98  }
99 
107  {
108  GEOS_UNUSED_VAR( stack );
109  return numNodes;
110  }
111 
113  virtual localIndex getMaxSupportPoints() const override
114  {
115  return maxSupportPoints;
116  }
117 
125  static void calcN( real64 const (&coords)[2],
126  real64 ( &N )[numNodes] );
127 
128 
137  static void calcN( localIndex const q,
138  real64 ( &N )[numNodes] );
139 
149  inline
150  static void calcN( localIndex const q,
151  StackVariables const & stack,
152  real64 ( &N )[numNodes] );
153 
161  static void calcBubbleN( real64 const (&pointCoord)[2],
162  real64 (& N)[1] )
163  {
164 
165  real64 const r = pointCoord[0];
166  real64 const s = pointCoord[1];
167 
168  N[0] = (1.0 - r - s) * r * s;
169 
170  }
171 
179  inline
180  static void calcBubbleN( localIndex const q,
181  real64 (& N)[1] )
182  {
183 
184  real64 const qCoords[2] = {quadratureParentCoords0( q ), quadratureParentCoords1( q ) };
185 
186  calcBubbleN( qCoords, N );
187  }
188 
198  real64 const (&X)[numNodes][3] );
199 
209  template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, bool UPPER >
211  inline
212  static void addGradGradStabilization( StackVariables const & stack,
213  real64 ( &matrix )
214  [maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT]
215  [maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT],
216  real64 const & scaleFactor );
217 
224  inline
225  static void getPermutation( int (& permutation)[numNodes] )
226  {
227  permutation[0] = 0;
228  permutation[1] = 1;
229  permutation[2] = 2;
230  }
231 
232 private:
234  constexpr static real64 parentArea = 0.5;
235 
237  //constexpr static real64 weight = parentArea / numQuadraturePoints;
238  constexpr static real64 weight = parentArea;
239 
246  inline
247  constexpr static real64 quadratureWeight( localIndex const q )
248  {
249 
250  if constexpr (numQuadraturePoints == 1)
251  {
252  constexpr real64 w[numQuadraturePoints] = { 1.0 };
253  return w[q];
254  }
255  else if constexpr (numQuadraturePoints == 4)
256  {
257  constexpr real64 w[numQuadraturePoints] = {-0.562500000000000,
258  0.520833333333333,
259  0.520833333333333,
260  0.520833333333333 };
261  return w[q];
262  }
263  else if constexpr (numQuadraturePoints == 6)
264  {
265  real64 const w[numQuadraturePoints] = { 0.166666666666666,
266  0.166666666666666,
267  0.166666666666666,
268  0.166666666666666,
269  0.166666666666666,
270  0.166666666666666 };
271  return w[q];
272  }
273 
274  }
275 
283  inline
284  constexpr static real64 quadratureParentCoords0( localIndex const q )
285  {
286 
287  if constexpr (numQuadraturePoints == 1)
288  {
289  constexpr real64 qCoords[numQuadraturePoints] = { 1.0/3.0 };
290  return qCoords[q];
291  }
292  else if constexpr (numQuadraturePoints == 4)
293  {
294  constexpr real64 qCoords[numQuadraturePoints] = { 0.333333333333333,
295  0.600000000000000,
296  0.200000000000000,
297  0.200000000000000 };
298  return qCoords[q];
299  }
300  else if constexpr (numQuadraturePoints == 6)
301  {
302  constexpr real64 qCoords[numQuadraturePoints] = { 0.659027622374092,
303  0.109039009072877,
304  0.231933368553031,
305  0.659027622374092,
306  0.109039009072877,
307  0.231933368553031 };
308  return qCoords[q];
309  }
310 
311  }
312 
320  inline
321  constexpr static real64 quadratureParentCoords1( localIndex const q )
322  {
323 
324  if constexpr (numQuadraturePoints == 1)
325  {
326  constexpr real64 qCoords[numQuadraturePoints] = { 1.0/3.0 };
327  return qCoords[q];
328  }
329  else if constexpr (numQuadraturePoints == 4)
330  {
331  constexpr real64 qCoords[numQuadraturePoints] = { 0.333333333333333,
332  0.200000000000000,
333  0.600000000000000,
334  0.200000000000000 };
335  return qCoords[q];
336  }
337  else if constexpr (numQuadraturePoints == 6)
338  {
339  constexpr real64 qCoords[numQuadraturePoints] = { 0.231933368553031,
340  0.659027622374092,
341  0.109039009072877,
342  0.109039009072877,
343  0.231933368553031,
344  0.659027622374092 };
345  return qCoords[q];
346  }
347 
348  }
349 
350 };
351 
353 
354 
355 template< typename NUM_Q_POINTS >
356 template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, bool UPPER >
358 inline
360 addGradGradStabilization( StackVariables const & stack,
361  real64 ( & matrix )
362  [maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT]
363  [maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT],
364  real64 const & scaleFactor )
365 {
366  GEOS_UNUSED_VAR( stack );
367  GEOS_UNUSED_VAR( matrix );
368  GEOS_UNUSED_VAR( scaleFactor );
369 }
370 
371 template< typename NUM_Q_POINTS >
373 inline
374 void
376 calcN( real64 const (&coords)[2],
377  real64 ( & N )[numNodes] )
378 {
379  real64 const r = coords[0];
380  real64 const s = coords[1];
381 
382  N[0] = 1.0 - r - s;
383  N[1] = r;
384  N[2] = s;
385 }
386 
387 template< typename NUM_Q_POINTS >
389 inline
390 void
392 calcN( localIndex const q,
393  real64 (& N)[numNodes] )
394 {
395 
396  real64 const qCoords[2] = {quadratureParentCoords0( q ), quadratureParentCoords1( q ) };
397 
398  calcN( qCoords, N );
399 
400 }
401 
402 template< typename NUM_Q_POINTS >
404 inline
406 calcN( localIndex const q,
407  StackVariables const & GEOS_UNUSED_PARAM( stack ),
408  real64 ( & N )[numNodes] )
409 {
410  return calcN( q, N );
411 }
412 
413 //*************************************************************************************************
414 
415 template< typename NUM_Q_POINTS >
417 inline
418 real64
421  real64 const (&X)[numNodes][3] )
422 {
423  //GEOS_UNUSED_VAR( q );
424  real64 n[3] = { ( X[1][1] - X[0][1] ) * ( X[2][2] - X[0][2] ) - ( X[2][1] - X[0][1] ) * ( X[1][2] - X[0][2] ),
425  ( X[2][0] - X[0][0] ) * ( X[1][2] - X[0][2] ) - ( X[1][0] - X[0][0] ) * ( X[2][2] - X[0][2] ),
426  ( X[1][0] - X[0][0] ) * ( X[2][1] - X[0][1] ) - ( X[2][0] - X[0][0] ) * ( X[1][1] - X[0][1] )};
427 
428  return sqrt( n[0] * n[0] + n[1] * n[1] + n[2] * n[2] ) * weight * quadratureWeight( q );
429 }
430 
432 
439 
440 
441 }
442 }
443 #endif //GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1TRIANGLEFACELAGRANGE1GAUSS_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.
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 calcBubbleN(localIndex const q, real64(&N)[1])
Calculate shape bubble functions values at a quadrature point.
constexpr static localIndex numNodes
The number of nodes/support points per element.
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.
constexpr static localIndex maxSupportPoints
The maximum 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 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 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3])
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.
static GEOS_HOST_DEVICE void getPermutation(int(&permutation)[numNodes])
Calculate the node permutation between the parent element and the geometric element....
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcN(real64 const (&coords)[2], real64(&N)[numNodes])
Calculate shape functions values at a single 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 getNumSupportPoints() const override
Virtual getter for the number of support points per element.
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...
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