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 >
107 template<
typename SUBREGION_TYPE >
222 template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS,
bool UPPER >
227 [MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT]
228 [MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT],
229 real64 const & scaleFactor )
237 for(
localIndex d = 0; d < NUMDOFSPERTRIALSUPPORTPOINT; ++d )
239 matrix[i*NUMDOFSPERTRIALSUPPORTPOINT + d][j*NUMDOFSPERTRIALSUPPORTPOINT + d] += contribution;
257 template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS >
261 real64 const ( &dofs )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
262 real64 ( & targetVector )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
263 real64 const scaleFactor )
265 for(
localIndex d = 0; d < NUMDOFSPERTRIALSUPPORTPOINT; ++d )
303 template<
typename SUBREGION_TYPE >
308 SUBREGION_TYPE
const & cellSubRegion,
313 meshData.
cellToNodeMap = cellSubRegion.nodeList().toViewConst();
314 meshData.
cellToFaceMap = cellSubRegion.faceList().toViewConst();
321 meshData.
cellCenters = cellSubRegion.getElementCenter();
322 meshData.
cellVolumes = cellSubRegion.getElementVolume();
332 template<
typename SUBREGION_TYPE >
344 computeProjectors< SUBREGION_TYPE >( cellIndex,
381 GEOS_ERROR(
"VEM functions have to be called with the StackVariables syntax" );
395 real64 (& samplingPointCoord)[3] )
412 GEOS_ERROR(
"VEM functions have to be called with the StackVariables syntax" );
431 GEOS_ERROR(
"VEM functions have to be called with the StackVariables syntax" );
453 GEOS_ERROR(
"No reference element map is defined for VEM classes" );
479 GEOS_ERROR(
"VEM functions have to be called with the StackVariables syntax" );
502 GEOS_ERROR(
"VEM functions have to be called with the StackVariables syntax" );
514 localIndex const (&faceToNodes)[MAXFACENODES],
515 localIndex const (&faceToEdges)[MAXFACENODES],
518 real64 const (&faceCenter)[3],
519 real64 const (&faceNormal)[3],
521 real64 const & invCellDiameter,
522 real64 const (&cellCenter)[3],
523 real64 ( &basisIntegrals )[MAXFACENODES],
524 real64 ( &threeDMonomialIntegrals )[3] );
535 template<
typename SUBREGION_TYPE >
539 computeProjectors(
localIndex const & cellIndex,
541 InputCellToNodeMap< SUBREGION_TYPE >
const & cellToNodeMap,
542 InputCellToFaceMap< SUBREGION_TYPE >
const & elementToFaceMap,
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]
558 template< localIndex DIMENSION,
typename POINT_COORDS_TYPE >
561 static real64 computeDiameter( POINT_COORDS_TYPE points,
565 for(
localIndex numPoint = 0; numPoint < numPoints; ++numPoint )
567 for(
localIndex numOthPoint = 0; numOthPoint < numPoint; ++numOthPoint )
569 real64 candidateDiameter = 0.0;
572 real64 coordDiff = points[numPoint][i] - points[numOthPoint][i];
573 candidateDiameter += coordDiff * coordDiff;
575 if( diameter < candidateDiameter )
577 diameter = candidateDiameter;
581 return LvArray::math::sqrt< real64 >( diameter );
584 template< localIndex DIMENSION,
typename POINT_COORDS_TYPE,
typename POINT_
SELECTION_TYPE >
587 static real64 computeDiameter( POINT_COORDS_TYPE points,
588 POINT_SELECTION_TYPE selectedPoints,
592 for(
localIndex numPoint = 0; numPoint < numSelectedPoints; ++numPoint )
594 for(
localIndex numOthPoint = 0; numOthPoint < numPoint; ++numOthPoint )
596 real64 candidateDiameter = 0.0;
599 real64 coordDiff = points[selectedPoints[numPoint]][i] -
600 points[selectedPoints[numOthPoint]][i];
601 candidateDiameter += coordDiff * coordDiff;
603 if( diameter < candidateDiameter )
605 diameter = candidateDiameter;
609 return LvArray::math::sqrt< real64 >( diameter );
615 #if !defined( GEOS_USE_HIP )
621 #if !defined( GEOS_USE_HIP )
638 #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.