20 #ifndef GEOS_VIRTUALELEMENT_CONFORMINGVIRTUALELEMENTORDER1_HPP_
21 #define GEOS_VIRTUALELEMENT_CONFORMINGVIRTUALELEMENTORDER1_HPP_
24 #include "codingUtilities/traits.hpp"
38 template< localIndex MAXCELLNODES, localIndex MAXFACENODES >
45 template<
typename SUBREGION_TYPE >
48 template<
typename SUBREGION_TYPE >
95 template<
typename SUBREGION_TYPE >
204 template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS,
bool UPPER >
209 [MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT]
210 [MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT],
211 real64 const & scaleFactor )
219 for(
localIndex d = 0; d < NUMDOFSPERTRIALSUPPORTPOINT; ++d )
221 matrix[i*NUMDOFSPERTRIALSUPPORTPOINT + d][j*NUMDOFSPERTRIALSUPPORTPOINT + d] += contribution;
239 template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS >
243 real64 const ( &dofs )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
244 real64 ( & targetVector )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
245 real64 const scaleFactor )
247 for(
localIndex d = 0; d < NUMDOFSPERTRIALSUPPORTPOINT; ++d )
285 template<
typename SUBREGION_TYPE >
290 SUBREGION_TYPE
const & cellSubRegion,
295 meshData.
cellToNodeMap = cellSubRegion.nodeList().toViewConst();
296 meshData.
cellToFaceMap = cellSubRegion.faceList().toViewConst();
303 meshData.
cellCenters = cellSubRegion.getElementCenter();
304 meshData.
cellVolumes = cellSubRegion.getElementVolume();
314 template<
typename SUBREGION_TYPE >
326 computeProjectors< SUBREGION_TYPE >( cellIndex,
363 GEOS_ERROR(
"VEM functions have to be called with the StackVariables syntax" );
377 real64 (& samplingPointCoord)[3] )
394 GEOS_ERROR(
"VEM functions have to be called with the StackVariables syntax" );
413 GEOS_ERROR(
"VEM functions have to be called with the StackVariables syntax" );
435 GEOS_ERROR(
"No reference element map is defined for VEM classes" );
461 GEOS_ERROR(
"VEM functions have to be called with the StackVariables syntax" );
484 GEOS_ERROR(
"VEM functions have to be called with the StackVariables syntax" );
496 localIndex const (&faceToNodes)[MAXFACENODES],
497 localIndex const (&faceToEdges)[MAXFACENODES],
500 real64 const (&faceCenter)[3],
501 real64 const (&faceNormal)[3],
503 real64 const & invCellDiameter,
504 real64 const (&cellCenter)[3],
505 real64 ( &basisIntegrals )[MAXFACENODES],
506 real64 ( &threeDMonomialIntegrals )[3] );
517 template<
typename SUBREGION_TYPE >
521 computeProjectors(
localIndex const & cellIndex,
523 InputCellToNodeMap< SUBREGION_TYPE >
const & cellToNodeMap,
524 InputCellToFaceMap< SUBREGION_TYPE >
const & elementToFaceMap,
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]
540 template< localIndex DIMENSION,
typename POINT_COORDS_TYPE >
543 static real64 computeDiameter( POINT_COORDS_TYPE points,
547 for(
localIndex numPoint = 0; numPoint < numPoints; ++numPoint )
549 for(
localIndex numOthPoint = 0; numOthPoint < numPoint; ++numOthPoint )
551 real64 candidateDiameter = 0.0;
554 real64 coordDiff = points[numPoint][i] - points[numOthPoint][i];
555 candidateDiameter += coordDiff * coordDiff;
557 if( diameter < candidateDiameter )
559 diameter = candidateDiameter;
563 return LvArray::math::sqrt< real64 >( diameter );
566 template< localIndex DIMENSION,
typename POINT_COORDS_TYPE,
typename POINT_
SELECTION_TYPE >
569 static real64 computeDiameter( POINT_COORDS_TYPE points,
570 POINT_SELECTION_TYPE selectedPoints,
574 for(
localIndex numPoint = 0; numPoint < numSelectedPoints; ++numPoint )
576 for(
localIndex numOthPoint = 0; numOthPoint < numPoint; ++numOthPoint )
578 real64 candidateDiameter = 0.0;
581 real64 coordDiff = points[selectedPoints[numPoint]][i] -
582 points[selectedPoints[numOthPoint]][i];
583 candidateDiameter += coordDiff * coordDiff;
585 if( diameter < candidateDiameter )
587 diameter = candidateDiameter;
591 return LvArray::math::sqrt< real64 >( diameter );
597 #if !defined( GEOS_USE_HIP )
603 #if !defined( GEOS_USE_HIP )
620 #if !defined( GEOS_USE_HIP )
#define GEOS_HOST_DEVICE
Marks a host-device function.
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
#define GEOS_FORCE_INLINE
Marks a function or lambda for inlining.
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
This class provides an interface to ObjectManagerBase in order to manage edge data.
NodeMapType & nodeList()
Get a node list.
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
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.
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.
LvArray::ArrayOfArraysView< T, INDEX_TYPE const, CONST_SIZES, LvArray::ChaiBuffer > ArrayOfArraysView
View of array of variable-sized arrays. See LvArray::ArrayOfArraysView for details.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Variables used to initialize the class.
Kernel variables allocated on the stack.