GEOS
FiniteElementBase.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 
21 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_FINITEELEMENTBASE_HPP_
22 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_FINITEELEMENTBASE_HPP_
23 
24 #include "common/DataTypes.hpp"
25 #include "common/GeosxMacros.hpp"
26 #include "finiteElement/PDEUtilities.hpp"
27 #include "LvArray/src/tensorOps.hpp"
28 #include "mesh/NodeManager.hpp"
29 #include "mesh/EdgeManager.hpp"
30 #include "mesh/FaceManager.hpp"
31 
32 namespace geos
33 {
34 namespace finiteElement
35 {
36 
38 template< int NUM_SUPPORT_POINTS,
39  int NUM_FACES,
40  int NUM_QUADRATURE_POINTS >
42 {
43 public:
44 
46  constexpr static localIndex numNodes = NUM_SUPPORT_POINTS;
47 
49  constexpr static localIndex numSupportPoints = NUM_SUPPORT_POINTS;
50 
53 
55  constexpr static localIndex numFaces = NUM_FACES;
56 
58  constexpr static localIndex numQuadraturePoints = NUM_QUADRATURE_POINTS;
59 
61  constexpr static int numSamplingPointsPerDirection = 10;
62 
65 
66 #ifdef __CUDACC__
67  #pragma diag_push
68  #pragma nv_diag_suppress 20012
69 #endif
88 #ifdef __CUDACC__
89  #pragma diag_pop
90 #endif
91 
100  {};
101 
106  template< typename SUBREGION_TYPE >
107  struct MeshData
108  {};
109 
116  {
117  return numQuadraturePoints;
118  }
119 
125  template< typename STACK_VARIABLES_TYPE >
127  static localIndex getNumQuadraturePoints( STACK_VARIABLES_TYPE const & stack )
128  {
129  GEOS_UNUSED_VAR( stack );
130  return numQuadraturePoints;
131  }
132 
139  {
140  return numNodes;
141  }
142 
148  template< typename STACK_VARIABLES_TYPE >
150  static localIndex getNumSupportPoints( STACK_VARIABLES_TYPE const & stack )
151  {
152  GEOS_UNUSED_VAR( stack );
153  return numNodes;
154  }
155 
162  {
163  return maxSupportPoints;
164  }
165 
174  static void getSamplingPointCoordInParentSpace( int const & linearIndex,
175  real64 (& samplingPointCoord)[3] )
176  {
177  GEOS_UNUSED_VAR( linearIndex, samplingPointCoord );
178  GEOS_ERROR( " Element type not supported." );
179  }
180 
181 
182 
191  template< typename SUBREGION_TYPE,
192  typename MESH_DATA_TYPE >
193  static void fillMeshData( NodeManager const & nodeManager,
194  EdgeManager const & edgeManager,
195  FaceManager const & faceManager,
196  SUBREGION_TYPE const & cellSubRegion,
197  MESH_DATA_TYPE & meshData )
198  {
199  GEOS_UNUSED_VAR( nodeManager,
200  edgeManager,
201  faceManager,
202  cellSubRegion,
203  meshData );
204  }
205 
216  template< typename LEAF,
217  typename SUBREGION_TYPE,
218  typename MESH_DATA_TYPE >
219  static void initialize( NodeManager const & nodeManager,
220  EdgeManager const & edgeManager,
221  FaceManager const & faceManager,
222  SUBREGION_TYPE const & cellSubRegion,
223  MESH_DATA_TYPE & meshData
224  )
225  {
226  LEAF::template fillMeshData< SUBREGION_TYPE >( nodeManager,
227  edgeManager,
228  faceManager,
229  cellSubRegion,
230  meshData );
231  }
232 
233 
240  template< typename MESH_DATA_TYPE,
241  typename STACK_VARIABLES_TYPE >
244  static void setupStack( localIndex const & cellIndex,
245  MESH_DATA_TYPE const & meshData,
246  STACK_VARIABLES_TYPE & stack )
247  {
248  GEOS_UNUSED_VAR( cellIndex,
249  meshData,
250  stack );
251 
252  }
253 
254 
263  template< typename LEAF,
264  typename MESH_DATA_TYPE >
266  void setup( localIndex const & cellIndex,
267  MESH_DATA_TYPE const & meshData,
268  typename LEAF::StackVariables & stack ) const
269  {
270  LEAF::setupStack( cellIndex, meshData, stack );
271  }
272 
273 
280  template< int N >
282  constexpr static PDEUtilities::FunctionSpace getFunctionSpace()
283  {
284  if constexpr ( N == 1 )
285  {
286  return PDEUtilities::FunctionSpace::H1;
287  }
288  else if constexpr ( N == 3 )
289  {
290  return PDEUtilities::FunctionSpace::H1vector;
291  }
292  else
293  {
294  static_assert( N == 1 || N == 3, "Unsupported number of components per support point" );
295  }
296  }
297 
298 
299 
310  template< localIndex NUMDOFSPERTRIALSUPPORTPOINT,
311  localIndex MAXSUPPORTPOINTS,
312  bool UPPER,
313  typename STACK_VARIABLES_TYPE >
316  static void addGradGradStabilization( STACK_VARIABLES_TYPE const & stack,
317  real64 ( & matrix )[MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT][MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT],
318  real64 const & scaleFactor )
319  {
320  GEOS_UNUSED_VAR( stack,
321  matrix,
322  scaleFactor );
323  }
324 
325 
326 
336  template< typename LEAF, localIndex NUMDOFSPERTRIALSUPPORTPOINT, bool UPPER = false >
338  void addGradGradStabilizationMatrix( typename LEAF::StackVariables const & stack,
339  real64 ( & matrix )[LEAF::maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT][LEAF::maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT],
340  real64 const scaleFactor = 1.0 ) const
341  {
342  LEAF::template addGradGradStabilization< NUMDOFSPERTRIALSUPPORTPOINT,
343  LEAF::maxSupportPoints,
344  UPPER >( stack,
345  matrix,
346  scaleFactor );
347  }
348 
362  template< localIndex NUMDOFSPERTRIALSUPPORTPOINT,
363  localIndex MAXSUPPORTPOINTS,
364  typename STACK_VARIABLES_TYPE >
367  static void addEvaluatedGradGradStabilization( STACK_VARIABLES_TYPE const & stack,
368  real64 const ( &dofs )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
369  real64 ( & targetVector )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
370  real64 const scaleFactor )
371  {
372  GEOS_UNUSED_VAR( stack );
373  GEOS_UNUSED_VAR( dofs );
374  GEOS_UNUSED_VAR( targetVector );
375  GEOS_UNUSED_VAR( scaleFactor );
376  }
377 
391  template< typename LEAF, localIndex NUMDOFSPERTRIALSUPPORTPOINT >
394  void
395  addEvaluatedGradGradStabilizationVector( typename LEAF::StackVariables const & stack,
396  real64 const ( &dofs )[LEAF::maxSupportPoints]
397  [NUMDOFSPERTRIALSUPPORTPOINT],
398  real64 ( & targetVector )[LEAF::maxSupportPoints]
399  [NUMDOFSPERTRIALSUPPORTPOINT],
400  real64 const scaleFactor = 1.0 ) const
401  {
402  LEAF::template
403  addEvaluatedGradGradStabilization< NUMDOFSPERTRIALSUPPORTPOINT, LEAF::maxSupportPoints >( stack,
404  dofs,
405  targetVector,
406  scaleFactor );
407  }
408 
409 
410 };
411 
412 
417 {
418 public:
419 
426  FiniteElementBase( localIndex const numSupportPoints,
427  localIndex const maxSupportPoints,
428  localIndex const numQuadraturePoints ):
429  m_numSupportPoints( numSupportPoints ),
430  m_maxSupportPoints( maxSupportPoints ),
431  m_numQuadraturePoints( numQuadraturePoints )
432  { }
433 
437  virtual ~FiniteElementBase() = default;
438 
444  {
445  return m_numQuadraturePoints;
446  }
447 
452  localIndex getNumSupportPoints() const { return m_numSupportPoints; };
453 
460  localIndex getMaxSupportPoints() const { return m_maxSupportPoints; };
461 
462 private:
463  localIndex const m_numSupportPoints;
464  localIndex const m_maxSupportPoints;
465  localIndex const m_numQuadraturePoints;
466 };
467 
468 } // namespace geos::finiteElement
469 } // namespace geos
470 
471 #endif //GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_FINITEELEMENTBASE_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_ERROR(...)
Raise a hard error and terminate the program.
Definition: Logger.hpp:204
This class provides an interface to ObjectManagerBase in order to manage edge data.
Definition: EdgeManager.hpp:43
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
Definition: FaceManager.hpp:44
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:46
Device-compatible (non virtual) Base class for all finite element formulations.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(STACK_VARIABLES_TYPE const &stack)
Get the number of quadrature points.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints()
Get the number of quadrature points.
constexpr static localIndex numQuadraturePoints
The number of quadrature points per element.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void addGradGradStabilization(STACK_VARIABLES_TYPE 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...
GEOS_HOST_DEVICE ~FiniteElementBase_impl()=default
Default destructor.
constexpr static int numSamplingPointsPerDirection
Number of sampling points.
GEOS_HOST_DEVICE FiniteElementBase_impl(FiniteElementBase_impl const &)=default
Default copy constructor.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints(STACK_VARIABLES_TYPE 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.
static void fillMeshData(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, SUBREGION_TYPE const &cellSubRegion, MESH_DATA_TYPE &meshData)
Method to fill a MeshData object.
GEOS_HOST_DEVICE void setup(localIndex const &cellIndex, MESH_DATA_TYPE const &meshData, typename LEAF::StackVariables &stack) const
Abstract setup method, possibly computing cell-dependent properties.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void addEvaluatedGradGradStabilization(STACK_VARIABLES_TYPE const &stack, real64 const (&dofs)[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT], real64(&targetVector)[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT], real64 const scaleFactor)
Empty method, here for compatibility with methods that require a stabilization of the grad-grad bilin...
GEOS_HOST_DEVICE FiniteElementBase_impl()=default
Default constructor.
constexpr static localIndex numSupportPoints
The number of support points per element.
GEOS_HOST_DEVICE FiniteElementBase_impl & operator=(FiniteElementBase_impl &&)=default
Default move assignment operator.
static GEOS_HOST_DEVICE localIndex getMaxSupportPoints()
Get the maximum number of support points.
constexpr static localIndex maxSupportPoints
The maximum number of support points per element.
GEOS_HOST_DEVICE void addGradGradStabilizationMatrix(typename LEAF::StackVariables const &stack, real64(&matrix)[LEAF::maxSupportPoints *NUMDOFSPERTRIALSUPPORTPOINT][LEAF::maxSupportPoints *NUMDOFSPERTRIALSUPPORTPOINT], real64 const scaleFactor=1.0) const
Add stabilization of grad-grad bilinear form to input matrix.
GEOS_HOST_DEVICE GEOS_FORCE_INLINE void addEvaluatedGradGradStabilizationVector(typename LEAF::StackVariables const &stack, real64 const (&dofs)[LEAF::maxSupportPoints][NUMDOFSPERTRIALSUPPORTPOINT], real64(&targetVector)[LEAF::maxSupportPoints][NUMDOFSPERTRIALSUPPORTPOINT], real64 const scaleFactor=1.0) const
Add a grad-grad stabilization operator evaluated at a provided vector of dofs to input vector.
constexpr static localIndex numNodes
The number of nodes per element.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints()
Get the number of support points.
constexpr static GEOS_HOST_DEVICE PDEUtilities::FunctionSpace getFunctionSpace()
Getter for the function space.
static void initialize(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, SUBREGION_TYPE const &cellSubRegion, MESH_DATA_TYPE &meshData)
Abstract initialization method.
constexpr static int numSamplingPoints
The number of sampling points per element.
GEOS_HOST_DEVICE FiniteElementBase_impl(FiniteElementBase_impl &&)=default
Default move constructor.
constexpr static localIndex numFaces
The number of faces per element.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void setupStack(localIndex const &cellIndex, MESH_DATA_TYPE const &meshData, STACK_VARIABLES_TYPE &stack)
Empty setup method.
GEOS_HOST_DEVICE FiniteElementBase_impl & operator=(FiniteElementBase_impl const &)=default
Default copy assignment operator.
Base class for FEM element implementations.
localIndex getNumSupportPoints() const
Getter for the number of support points per element.
virtual ~FiniteElementBase()=default
Destructor.
localIndex getNumQuadraturePoints() const
Getter for the number of quadrature points per element.
localIndex getMaxSupportPoints() const
Get the maximum number of support points for this element.
FiniteElementBase(localIndex const numSupportPoints, localIndex const maxSupportPoints, localIndex const numQuadraturePoints)
Default constructor.
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