GEOS
ConformingVirtualElementOrder1.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_VIRTUALELEMENT_CONFORMINGVIRTUALELEMENTORDER1_HPP_
21 #define GEOS_VIRTUALELEMENT_CONFORMINGVIRTUALELEMENTORDER1_HPP_
22 
23 #include "FiniteElementBase.hpp"
24 #include "codingUtilities/traits.hpp"
25 
26 namespace geos
27 {
28 namespace finiteElement
29 {
38 template< localIndex MAXCELLNODES, localIndex MAXFACENODES >
40 {
41 public:
45  template< typename SUBREGION_TYPE >
46  using InputCellToNodeMap = traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType >;
48  template< typename SUBREGION_TYPE >
49  using InputCellToFaceMap = traits::ViewTypeConst< typename SUBREGION_TYPE::FaceMapType >;
56 
58  static constexpr localIndex maxSupportPoints = MAXCELLNODES;
60  static constexpr localIndex numNodes = MAXCELLNODES;
62  static constexpr localIndex numQuadraturePoints = 1;
63 
73  {
74 
84  real64 stabilizationMatrix[MAXCELLNODES][MAXCELLNODES];
88  };
89 
95  template< typename SUBREGION_TYPE >
96  struct MeshData : public FiniteElementBase::MeshData< SUBREGION_TYPE >
97  {
120  };
121 
123  inline
125  {
126  return numQuadraturePoints;
127  }
128 
130  inline
131  virtual localIndex getMaxSupportPoints() const override
132  {
133  return maxSupportPoints;
134  }
135 
142  inline
144  {
145  return stack.numSupportPoints;
146  }
147 
159  inline
160  static real64 calcGradN( localIndex const q,
161  real64 const (&X)[maxSupportPoints][3],
162  StackVariables const & stack,
163  real64 ( & gradN )[maxSupportPoints][3] )
164  {
165  for( localIndex i = 0; i < stack.numSupportPoints; ++i )
166  {
167  gradN[i][0] = stack.basisDerivativesIntegralMean[i][0];
168  gradN[i][1] = stack.basisDerivativesIntegralMean[i][1];
169  gradN[i][2] = stack.basisDerivativesIntegralMean[i][2];
170  }
171  return transformedQuadratureWeight( q, X, stack );
172  }
173 
183  inline
184  static void calcN( localIndex const q,
185  StackVariables const & stack,
186  real64 ( & N )[maxSupportPoints] )
187  {
188  GEOS_UNUSED_VAR( q );
189  for( localIndex i = 0; i < stack.numSupportPoints; ++i )
190  {
191  N[i] = stack.basisFunctionsIntegralMean[i];
192  }
193  }
194 
204  template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS, bool UPPER >
206  inline
207  static void addGradGradStabilization( StackVariables const & stack,
208  real64 ( & matrix )
209  [MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT]
210  [MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT],
211  real64 const & scaleFactor )
212  {
213  for( localIndex i = 0; i < stack.numSupportPoints; ++i )
214  {
215  localIndex startCol = (UPPER) ? i : 0;
216  for( localIndex j = startCol; j < stack.numSupportPoints; ++j )
217  {
218  real64 const contribution = scaleFactor * stack.stabilizationMatrix[i][j];
219  for( localIndex d = 0; d < NUMDOFSPERTRIALSUPPORTPOINT; ++d )
220  {
221  matrix[i*NUMDOFSPERTRIALSUPPORTPOINT + d][j*NUMDOFSPERTRIALSUPPORTPOINT + d] += contribution;
222  }
223  }
224  }
225  }
226 
239  template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS >
241  inline
243  real64 const ( &dofs )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
244  real64 ( & targetVector )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
245  real64 const scaleFactor )
246  {
247  for( localIndex d = 0; d < NUMDOFSPERTRIALSUPPORTPOINT; ++d )
248  {
249  for( localIndex i = 0; i < stack.numSupportPoints; ++i )
250  {
251  for( localIndex j = 0; j < stack.numSupportPoints; ++j )
252  {
253  targetVector[i][d] += scaleFactor * stack.stabilizationMatrix[i][j] * dofs[j][d];
254  }
255  }
256  }
257  }
258 
268  inline
270  real64 const ( &X )[maxSupportPoints][3],
271  StackVariables const & stack )
272  {
273  GEOS_UNUSED_VAR( q, X );
274  return stack.quadratureWeight;
275  }
276 
285  template< typename SUBREGION_TYPE >
286  inline
287  static void fillMeshData( NodeManager const & nodeManager,
288  EdgeManager const & edgeManager,
289  FaceManager const & faceManager,
290  SUBREGION_TYPE const & cellSubRegion,
291  MeshData< SUBREGION_TYPE > & meshData
292  )
293  {
294  meshData.nodesCoords = nodeManager.referencePosition();
295  meshData.cellToNodeMap = cellSubRegion.nodeList().toViewConst();
296  meshData.cellToFaceMap = cellSubRegion.faceList().toViewConst();
297  meshData.faceToNodeMap = faceManager.nodeList().toViewConst();
298  meshData.faceToEdgeMap = faceManager.edgeList().toViewConst();
299  meshData.edgeToNodeMap = edgeManager.nodeList().toViewConst();
300  meshData.faceCenters = faceManager.faceCenter();
301  meshData.faceNormals = faceManager.faceNormal();
302  meshData.faceAreas = faceManager.faceArea();
303  meshData.cellCenters = cellSubRegion.getElementCenter();
304  meshData.cellVolumes = cellSubRegion.getElementVolume();
305  }
306 
314  template< typename SUBREGION_TYPE >
316  inline
317  static void setupStack( localIndex const & cellIndex,
318  MeshData< SUBREGION_TYPE > const & meshData,
319  StackVariables & stack )
320  {
321  real64 const cellCenter[3] { meshData.cellCenters( cellIndex, 0 ),
322  meshData.cellCenters( cellIndex, 1 ),
323  meshData.cellCenters( cellIndex, 2 ) };
324  real64 const cellVolume = meshData.cellVolumes( cellIndex );
325 
326  computeProjectors< SUBREGION_TYPE >( cellIndex,
327  meshData.nodesCoords,
328  meshData.cellToNodeMap,
329  meshData.cellToFaceMap,
330  meshData.faceToNodeMap,
331  meshData.faceToEdgeMap,
332  meshData.edgeToNodeMap,
333  meshData.faceCenters,
334  meshData.faceNormals,
335  meshData.faceAreas,
336  cellCenter,
337  cellVolume,
338  stack.numSupportPoints,
339  stack.quadratureWeight,
341  stack.stabilizationMatrix,
343  }
344 
360  inline
362  {
363  GEOS_ERROR( "VEM functions have to be called with the StackVariables syntax" );
364  return 0;
365  }
366 
367 
376  static void getSamplingPointCoordInParentSpace( int const & linearIndex,
377  real64 (& samplingPointCoord)[3] )
378  {
379  GEOS_UNUSED_VAR( linearIndex, samplingPointCoord );
380  GEOS_ERROR( "Element type not supported." );
381  }
382 
390  inline
391  static void calcN( localIndex const q,
392  real64 ( & N )[maxSupportPoints] )
393  {
394  GEOS_ERROR( "VEM functions have to be called with the StackVariables syntax" );
395  GEOS_UNUSED_VAR( q );
396  for( localIndex i = 0; i < maxSupportPoints; ++i )
397  {
398  N[i] = 0.0;
399  }
400  }
401 
409  inline
410  static void calcN( real64 const (&coords)[3],
411  real64 (& N)[maxSupportPoints] )
412  {
413  GEOS_ERROR( "VEM functions have to be called with the StackVariables syntax" );
414  GEOS_UNUSED_VAR( coords );
415  for( localIndex i = 0; i < maxSupportPoints; ++i )
416  {
417  N[i] = 0.0;
418  }
419  }
420 
430  inline
431  static real64 invJacobianTransformation( int const q,
432  real64 const ( &X )[numNodes][3],
433  real64 ( & J )[3][3] )
434  {
435  GEOS_ERROR( "No reference element map is defined for VEM classes" );
436  GEOS_UNUSED_VAR( q, X );
437  for( localIndex i = 0; i < 3; ++i )
438  {
439  for( localIndex j = 0; j < 3; ++j )
440  {
441  J[i][j] = 0.0;
442  }
443  }
444  return 0.0;
445  }
446 
456  inline
457  static real64 calcGradN( localIndex const q,
458  real64 const ( &X )[maxSupportPoints][3],
459  real64 ( & gradN )[maxSupportPoints][3] )
460  {
461  GEOS_ERROR( "VEM functions have to be called with the StackVariables syntax" );
462  GEOS_UNUSED_VAR( q, X );
463  for( localIndex i = 0; i < maxSupportPoints; ++i )
464  {
465  for( localIndex j = 0; j < 3; ++j )
466  {
467  gradN[i][j] = 0.0;
468  }
469  }
470  return 0.0;
471  }
472 
482  real64 const ( &X )[maxSupportPoints][3] ) const
483  {
484  GEOS_ERROR( "VEM functions have to be called with the StackVariables syntax" );
485  GEOS_UNUSED_VAR( q, X );
486  return 0.0;
487  }
488 
491 private:
492 
494  static void
495  computeFaceIntegrals( InputNodeCoords const & nodesCoords,
496  localIndex const (&faceToNodes)[MAXFACENODES],
497  localIndex const (&faceToEdges)[MAXFACENODES],
498  localIndex const & numFaceVertices,
499  real64 const & faceArea,
500  real64 const (&faceCenter)[3],
501  real64 const (&faceNormal)[3],
502  InputEdgeToNodeMap const & edgeToNodes,
503  real64 const & invCellDiameter,
504  real64 const (&cellCenter)[3],
505  real64 ( &basisIntegrals )[MAXFACENODES],
506  real64 ( &threeDMonomialIntegrals )[3] );
507 
517  template< typename SUBREGION_TYPE >
519  inline
520  static void
521  computeProjectors( localIndex const & cellIndex,
522  InputNodeCoords const & nodesCoords,
523  InputCellToNodeMap< SUBREGION_TYPE > const & cellToNodeMap,
524  InputCellToFaceMap< SUBREGION_TYPE > const & elementToFaceMap,
525  InputFaceToNodeMap const & faceToNodeMap,
526  InputFaceToEdgeMap const & faceToEdgeMap,
527  InputEdgeToNodeMap const & edgeToNodeMap,
528  arrayView2d< real64 const > const faceCenters,
529  arrayView2d< real64 const > const faceNormals,
530  arrayView1d< real64 const > const faceAreas,
531  real64 const (&cellCenter)[3],
532  real64 const & cellVolume,
534  real64 & quadratureWeight,
535  real64 ( &basisFunctionsIntegralMean )[MAXCELLNODES],
536  real64 ( &stabilizationMatrix )[MAXCELLNODES][MAXCELLNODES],
537  real64 ( &basisDerivativesIntegralMean )[MAXCELLNODES][3]
538  );
539 
540  template< localIndex DIMENSION, typename POINT_COORDS_TYPE >
542  inline
543  static real64 computeDiameter( POINT_COORDS_TYPE points,
544  localIndex const & numPoints )
545  {
546  real64 diameter = 0;
547  for( localIndex numPoint = 0; numPoint < numPoints; ++numPoint )
548  {
549  for( localIndex numOthPoint = 0; numOthPoint < numPoint; ++numOthPoint )
550  {
551  real64 candidateDiameter = 0.0;
552  for( localIndex i = 0; i < DIMENSION; ++i )
553  {
554  real64 coordDiff = points[numPoint][i] - points[numOthPoint][i];
555  candidateDiameter += coordDiff * coordDiff;
556  }
557  if( diameter < candidateDiameter )
558  {
559  diameter = candidateDiameter;
560  }
561  }
562  }
563  return LvArray::math::sqrt< real64 >( diameter );
564  }
565 
566  template< localIndex DIMENSION, typename POINT_COORDS_TYPE, typename POINT_SELECTION_TYPE >
568  inline
569  static real64 computeDiameter( POINT_COORDS_TYPE points,
570  POINT_SELECTION_TYPE selectedPoints,
571  localIndex const & numSelectedPoints )
572  {
573  real64 diameter = 0;
574  for( localIndex numPoint = 0; numPoint < numSelectedPoints; ++numPoint )
575  {
576  for( localIndex numOthPoint = 0; numOthPoint < numPoint; ++numOthPoint )
577  {
578  real64 candidateDiameter = 0.0;
579  for( localIndex i = 0; i < DIMENSION; ++i )
580  {
581  real64 coordDiff = points[selectedPoints[numPoint]][i] -
582  points[selectedPoints[numOthPoint]][i];
583  candidateDiameter += coordDiff * coordDiff;
584  }
585  if( diameter < candidateDiameter )
586  {
587  diameter = candidateDiameter;
588  }
589  }
590  }
591  return LvArray::math::sqrt< real64 >( diameter );
592  }
593 };
594 
597 #if !defined( GEOS_USE_HIP )
600 #endif
603 #if !defined( GEOS_USE_HIP )
606 #endif
620 #if !defined( GEOS_USE_HIP )
622 #endif
623 }
624 }
625 
627 
628 #endif // GEOS_VIRTUALELEMENT_CONFORMINGVIRTUALELEMENTORDER1_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(msg)
Raise a hard error and terminate the program.
Definition: Logger.hpp:157
This class provides an interface to ObjectManagerBase in order to manage edge data.
Definition: EdgeManager.hpp:43
NodeMapType & nodeList()
Get a node list.
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
Definition: FaceManager.hpp:44
array2d< real64 > & faceNormal()
Get a mutable accessor to an array containing all the face normals.
array1d< real64 > & faceArea()
Get a mutable accessor to an array containing all the face area.
array2d< real64 > & faceCenter()
Get a mutable accessor to an array containing all the face center.
NodeMapType & nodeList()
Get a mutable accessor to a map containing the list of each nodes for each faces.
EdgeMapType & edgeList()
Get a mutable accessor to a map containing the list of each edges for each faces.
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:46
GEOS_HOST_DEVICE localIndex getNumQuadraturePoints() const override
Virtual getter for the number of quadrature points per element.
ArrayOfArraysView< localIndex const > InputFaceToNodeMap
Type of MeshData::faceToNodeMap.
static constexpr localIndex maxSupportPoints
The maximum number of support points per element.
static GEOS_HOST_DEVICE void addEvaluatedGradGradStabilization(StackVariables const &stack, real64 const (&dofs)[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT], real64(&targetVector)[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT], real64 const scaleFactor)
Adds a grad-grad stabilization evaluated at dofs to targetVector.
arrayView2d< localIndex const > InputEdgeToNodeMap
Type of MeshData::edgeToNodeMap.
static GEOS_HOST_DEVICE real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[maxSupportPoints][3], StackVariables const &stack)
Calculate the integration weights for a quadrature point.
static constexpr localIndex numQuadraturePoints
The number of quadrature points per element.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints(StackVariables const &stack)
Get the number of support points.
ArrayOfArraysView< localIndex const > InputFaceToEdgeMap
Type of MeshData::faceToEdgeMap.
traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType > InputCellToNodeMap
Type of MeshData::cellToNodeMap.
static void fillMeshData(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, SUBREGION_TYPE const &cellSubRegion, MeshData< SUBREGION_TYPE > &meshData)
Method to fill a MeshData object.
static GEOS_HOST_DEVICE real64 calcGradN(localIndex const q, real64 const (&X)[maxSupportPoints][3], StackVariables const &stack, real64(&gradN)[maxSupportPoints][3])
Calculate the shape functions projected derivatives wrt the physical coordinates.
static GEOS_HOST_DEVICE void addGradGradStabilization(StackVariables const &stack, real64(&matrix)[MAXSUPPORTPOINTS *NUMDOFSPERTRIALSUPPORTPOINT][MAXSUPPORTPOINTS *NUMDOFSPERTRIALSUPPORTPOINT], real64 const &scaleFactor)
Adds a grad-grad stabilization to matrix.
static constexpr localIndex numNodes
Static property kept for consistency with other finite element classes.
traits::ViewTypeConst< typename SUBREGION_TYPE::FaceMapType > InputCellToFaceMap
Type of MeshData::cellToFaceMap.
virtual GEOS_HOST_DEVICE localIndex getMaxSupportPoints() const override
Get the maximum number of support points for this element.
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > InputNodeCoords
Type of MeshData::nodesCoords.
static GEOS_HOST_DEVICE void calcN(localIndex const q, StackVariables const &stack, real64(&N)[maxSupportPoints])
Calculate shape functions values for each support point at a quadrature point.
static GEOS_HOST_DEVICE void setupStack(localIndex const &cellIndex, MeshData< SUBREGION_TYPE > const &meshData, StackVariables &stack)
Setup method.
Base class for FEM element implementations.
GEOS_HOST_DEVICE localIndex numSupportPoints(typename LEAF::StackVariables const &stack) const
Getter for the number of support points per element.
array2d< real64, nodes::REFERENCE_POSITION_PERM > & referencePosition()
Get the mutable reference position array. This table will contain all the node coordinates.
static GEOS_HOST_DEVICE real64 invJacobianTransformation(int const q, real64 const (&X)[numNodes][3], real64(&J)[3][3])
This function returns an error, since there is no reference element defined in the VEM context....
static GEOS_HOST_DEVICE real64 calcGradN(localIndex const q, real64 const (&X)[maxSupportPoints][3], real64(&gradN)[maxSupportPoints][3])
This function returns an error, since to get projection of gradients with VEM you have to use the Sta...
GEOS_HOST_DEVICE real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[maxSupportPoints][3]) const
This function returns an error, since to get the quadrature weight with VEM you have to use the Stack...
GEOS_HOST_DEVICE localIndex getNumSupportPoints() const override
This function returns an error, since to get the number of support points with VEM you have to use th...
static GEOS_HOST_DEVICE void calcN(real64 const (&coords)[3], real64(&N)[maxSupportPoints])
This function returns an error, since to get projection of gradients with VEM you have to use the Sta...
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 void calcN(localIndex const q, real64(&N)[maxSupportPoints])
This function returns an error, since to get projection of basis functions with VEM you have to use t...
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
LvArray::ArrayOfArraysView< T, INDEX_TYPE const, CONST_SIZES, LvArray::ChaiBuffer > ArrayOfArraysView
View of array of variable-sized arrays. See LvArray::ArrayOfArraysView for details.
Definition: DataTypes.hpp:286
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
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:196
InputFaceToNodeMap faceToNodeMap
View to the face-to-node map in the sub-region.
InputEdgeToNodeMap edgeToNodeMap
View to the edge-to-node map in the sub-region.
arrayView1d< real64 const > faceAreas
View to the array of face areas.
arrayView2d< real64 const > faceNormals
View to the array of face normals.
InputCellToNodeMap< SUBREGION_TYPE > cellToNodeMap
View to the cell-to-node map in the sub-region.
InputFaceToEdgeMap faceToEdgeMap
View to the face-to-edge map in the sub-region.
arrayView2d< real64 const > faceCenters
View to the array of face centers.
arrayView1d< real64 const > cellVolumes
View to the array of cell volumes.
InputCellToFaceMap< SUBREGION_TYPE > cellToFaceMap
View to the cell-to-face map in the sub-region.
InputNodeCoords nodesCoords
View to the array containing nodes coordinates.
arrayView2d< real64 const > cellCenters
View to the array of cell centers.
Variables used to initialize the class.