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 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_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 
65 
67  virtual ~ConformingVirtualElementOrder1() override
68  {}
69 
79  {
85  {}
86 
96  real64 stabilizationMatrix[MAXCELLNODES][MAXCELLNODES];
100  };
101 
107  template< typename SUBREGION_TYPE >
108  struct MeshData : public FiniteElementBase::MeshData< SUBREGION_TYPE >
109  {
114  {}
115 
138  };
139 
141  inline
143  {
144  return numQuadraturePoints;
145  }
146 
148  inline
149  virtual localIndex getMaxSupportPoints() const override
150  {
151  return maxSupportPoints;
152  }
153 
160  inline
162  {
163  return stack.numSupportPoints;
164  }
165 
177  inline
178  static real64 calcGradN( localIndex const q,
179  real64 const (&X)[maxSupportPoints][3],
180  StackVariables const & stack,
181  real64 ( & gradN )[maxSupportPoints][3] )
182  {
183  for( localIndex i = 0; i < stack.numSupportPoints; ++i )
184  {
185  gradN[i][0] = stack.basisDerivativesIntegralMean[i][0];
186  gradN[i][1] = stack.basisDerivativesIntegralMean[i][1];
187  gradN[i][2] = stack.basisDerivativesIntegralMean[i][2];
188  }
189  return transformedQuadratureWeight( q, X, stack );
190  }
191 
201  inline
202  static void calcN( localIndex const q,
203  StackVariables const & stack,
204  real64 ( & N )[maxSupportPoints] )
205  {
206  GEOS_UNUSED_VAR( q );
207  for( localIndex i = 0; i < stack.numSupportPoints; ++i )
208  {
209  N[i] = stack.basisFunctionsIntegralMean[i];
210  }
211  }
212 
222  template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS, bool UPPER >
224  inline
225  static void addGradGradStabilization( StackVariables const & stack,
226  real64 ( & matrix )
227  [MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT]
228  [MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT],
229  real64 const & scaleFactor )
230  {
231  for( localIndex i = 0; i < stack.numSupportPoints; ++i )
232  {
233  localIndex startCol = (UPPER) ? i : 0;
234  for( localIndex j = startCol; j < stack.numSupportPoints; ++j )
235  {
236  real64 const contribution = scaleFactor * stack.stabilizationMatrix[i][j];
237  for( localIndex d = 0; d < NUMDOFSPERTRIALSUPPORTPOINT; ++d )
238  {
239  matrix[i*NUMDOFSPERTRIALSUPPORTPOINT + d][j*NUMDOFSPERTRIALSUPPORTPOINT + d] += contribution;
240  }
241  }
242  }
243  }
244 
257  template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS >
259  inline
261  real64 const ( &dofs )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
262  real64 ( & targetVector )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
263  real64 const scaleFactor )
264  {
265  for( localIndex d = 0; d < NUMDOFSPERTRIALSUPPORTPOINT; ++d )
266  {
267  for( localIndex i = 0; i < stack.numSupportPoints; ++i )
268  {
269  for( localIndex j = 0; j < stack.numSupportPoints; ++j )
270  {
271  targetVector[i][d] += scaleFactor * stack.stabilizationMatrix[i][j] * dofs[j][d];
272  }
273  }
274  }
275  }
276 
286  inline
288  real64 const ( &X )[maxSupportPoints][3],
289  StackVariables const & stack )
290  {
291  GEOS_UNUSED_VAR( q, X );
292  return stack.quadratureWeight;
293  }
294 
303  template< typename SUBREGION_TYPE >
304  inline
305  static void fillMeshData( NodeManager const & nodeManager,
306  EdgeManager const & edgeManager,
307  FaceManager const & faceManager,
308  SUBREGION_TYPE const & cellSubRegion,
309  MeshData< SUBREGION_TYPE > & meshData
310  )
311  {
312  meshData.nodesCoords = nodeManager.referencePosition();
313  meshData.cellToNodeMap = cellSubRegion.nodeList().toViewConst();
314  meshData.cellToFaceMap = cellSubRegion.faceList().toViewConst();
315  meshData.faceToNodeMap = faceManager.nodeList().toViewConst();
316  meshData.faceToEdgeMap = faceManager.edgeList().toViewConst();
317  meshData.edgeToNodeMap = edgeManager.nodeList().toViewConst();
318  meshData.faceCenters = faceManager.faceCenter();
319  meshData.faceNormals = faceManager.faceNormal();
320  meshData.faceAreas = faceManager.faceArea();
321  meshData.cellCenters = cellSubRegion.getElementCenter();
322  meshData.cellVolumes = cellSubRegion.getElementVolume();
323  }
324 
332  template< typename SUBREGION_TYPE >
334  inline
335  static void setupStack( localIndex const & cellIndex,
336  MeshData< SUBREGION_TYPE > const & meshData,
337  StackVariables & stack )
338  {
339  real64 const cellCenter[3] { meshData.cellCenters( cellIndex, 0 ),
340  meshData.cellCenters( cellIndex, 1 ),
341  meshData.cellCenters( cellIndex, 2 ) };
342  real64 const cellVolume = meshData.cellVolumes( cellIndex );
343 
344  computeProjectors< SUBREGION_TYPE >( cellIndex,
345  meshData.nodesCoords,
346  meshData.cellToNodeMap,
347  meshData.cellToFaceMap,
348  meshData.faceToNodeMap,
349  meshData.faceToEdgeMap,
350  meshData.edgeToNodeMap,
351  meshData.faceCenters,
352  meshData.faceNormals,
353  meshData.faceAreas,
354  cellCenter,
355  cellVolume,
356  stack.numSupportPoints,
357  stack.quadratureWeight,
359  stack.stabilizationMatrix,
361  }
362 
378  inline
380  {
381  GEOS_ERROR( "VEM functions have to be called with the StackVariables syntax" );
382  return 0;
383  }
384 
385 
394  static void getSamplingPointCoordInParentSpace( int const & linearIndex,
395  real64 (& samplingPointCoord)[3] )
396  {
397  GEOS_UNUSED_VAR( linearIndex, samplingPointCoord );
398  GEOS_ERROR( "Element type not supported." );
399  }
400 
408  inline
409  static void calcN( localIndex const q,
410  real64 ( & N )[maxSupportPoints] )
411  {
412  GEOS_ERROR( "VEM functions have to be called with the StackVariables syntax" );
413  GEOS_UNUSED_VAR( q );
414  for( localIndex i = 0; i < maxSupportPoints; ++i )
415  {
416  N[i] = 0.0;
417  }
418  }
419 
427  inline
428  static void calcN( real64 const (&coords)[3],
429  real64 (& N)[maxSupportPoints] )
430  {
431  GEOS_ERROR( "VEM functions have to be called with the StackVariables syntax" );
432  GEOS_UNUSED_VAR( coords );
433  for( localIndex i = 0; i < maxSupportPoints; ++i )
434  {
435  N[i] = 0.0;
436  }
437  }
438 
448  inline
449  static real64 invJacobianTransformation( int const q,
450  real64 const ( &X )[numNodes][3],
451  real64 ( & J )[3][3] )
452  {
453  GEOS_ERROR( "No reference element map is defined for VEM classes" );
454  GEOS_UNUSED_VAR( q, X );
455  for( localIndex i = 0; i < 3; ++i )
456  {
457  for( localIndex j = 0; j < 3; ++j )
458  {
459  J[i][j] = 0.0;
460  }
461  }
462  return 0.0;
463  }
464 
474  inline
475  static real64 calcGradN( localIndex const q,
476  real64 const ( &X )[maxSupportPoints][3],
477  real64 ( & gradN )[maxSupportPoints][3] )
478  {
479  GEOS_ERROR( "VEM functions have to be called with the StackVariables syntax" );
480  GEOS_UNUSED_VAR( q, X );
481  for( localIndex i = 0; i < maxSupportPoints; ++i )
482  {
483  for( localIndex j = 0; j < 3; ++j )
484  {
485  gradN[i][j] = 0.0;
486  }
487  }
488  return 0.0;
489  }
490 
500  real64 const ( &X )[maxSupportPoints][3] ) const
501  {
502  GEOS_ERROR( "VEM functions have to be called with the StackVariables syntax" );
503  GEOS_UNUSED_VAR( q, X );
504  return 0.0;
505  }
506 
509 private:
510 
512  static void
513  computeFaceIntegrals( InputNodeCoords const & nodesCoords,
514  localIndex const (&faceToNodes)[MAXFACENODES],
515  localIndex const (&faceToEdges)[MAXFACENODES],
516  localIndex const & numFaceVertices,
517  real64 const & faceArea,
518  real64 const (&faceCenter)[3],
519  real64 const (&faceNormal)[3],
520  InputEdgeToNodeMap const & edgeToNodes,
521  real64 const & invCellDiameter,
522  real64 const (&cellCenter)[3],
523  real64 ( &basisIntegrals )[MAXFACENODES],
524  real64 ( &threeDMonomialIntegrals )[3] );
525 
535  template< typename SUBREGION_TYPE >
537  inline
538  static void
539  computeProjectors( localIndex const & cellIndex,
540  InputNodeCoords const & nodesCoords,
541  InputCellToNodeMap< SUBREGION_TYPE > const & cellToNodeMap,
542  InputCellToFaceMap< SUBREGION_TYPE > const & elementToFaceMap,
543  InputFaceToNodeMap const & faceToNodeMap,
544  InputFaceToEdgeMap const & faceToEdgeMap,
545  InputEdgeToNodeMap const & edgeToNodeMap,
546  arrayView2d< real64 const > const faceCenters,
547  arrayView2d< real64 const > const faceNormals,
548  arrayView1d< real64 const > const faceAreas,
549  real64 const (&cellCenter)[3],
550  real64 const & cellVolume,
552  real64 & quadratureWeight,
553  real64 ( &basisFunctionsIntegralMean )[MAXCELLNODES],
554  real64 ( &stabilizationMatrix )[MAXCELLNODES][MAXCELLNODES],
555  real64 ( &basisDerivativesIntegralMean )[MAXCELLNODES][3]
556  );
557 
558  template< localIndex DIMENSION, typename POINT_COORDS_TYPE >
560  inline
561  static real64 computeDiameter( POINT_COORDS_TYPE points,
562  localIndex const & numPoints )
563  {
564  real64 diameter = 0;
565  for( localIndex numPoint = 0; numPoint < numPoints; ++numPoint )
566  {
567  for( localIndex numOthPoint = 0; numOthPoint < numPoint; ++numOthPoint )
568  {
569  real64 candidateDiameter = 0.0;
570  for( localIndex i = 0; i < DIMENSION; ++i )
571  {
572  real64 coordDiff = points[numPoint][i] - points[numOthPoint][i];
573  candidateDiameter += coordDiff * coordDiff;
574  }
575  if( diameter < candidateDiameter )
576  {
577  diameter = candidateDiameter;
578  }
579  }
580  }
581  return LvArray::math::sqrt< real64 >( diameter );
582  }
583 
584  template< localIndex DIMENSION, typename POINT_COORDS_TYPE, typename POINT_SELECTION_TYPE >
586  inline
587  static real64 computeDiameter( POINT_COORDS_TYPE points,
588  POINT_SELECTION_TYPE selectedPoints,
589  localIndex const & numSelectedPoints )
590  {
591  real64 diameter = 0;
592  for( localIndex numPoint = 0; numPoint < numSelectedPoints; ++numPoint )
593  {
594  for( localIndex numOthPoint = 0; numOthPoint < numPoint; ++numOthPoint )
595  {
596  real64 candidateDiameter = 0.0;
597  for( localIndex i = 0; i < DIMENSION; ++i )
598  {
599  real64 coordDiff = points[selectedPoints[numPoint]][i] -
600  points[selectedPoints[numOthPoint]][i];
601  candidateDiameter += coordDiff * coordDiff;
602  }
603  if( diameter < candidateDiameter )
604  {
605  diameter = candidateDiameter;
606  }
607  }
608  }
609  return LvArray::math::sqrt< real64 >( diameter );
610  }
611 };
612 
615 #if !defined( GEOS_USE_HIP )
618 #endif
621 #if !defined( GEOS_USE_HIP )
624 #endif
638 #if !defined( GEOS_USE_HIP )
640 #endif
641 }
642 }
643 
645 
646 #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.