20 #ifndef GEOS_MESH_UTILITIES_CICOMPUTATIONKERNEL_HPP_ 
   21 #define GEOS_MESH_UTILITIES_CICOMPUTATIONKERNEL_HPP_ 
   25 #include "finiteElement/FiniteElementDispatch.hpp" 
   26 #include "common/GEOS_RAJA_Interface.hpp" 
   28 #include "mesh/CellElementSubRegion.hpp" 
   40 template< 
typename FE_TYPE >
 
   56     m_finiteElementSpace( finiteElementSpace ),
 
   57     m_elementType( elementSubRegion.getElementType() ),
 
   58     m_X( nodeManager.referencePosition() ),
 
   59     m_elemsToNodes( elementSubRegion.nodeList().toViewConst() ),
 
   60     m_fracturedElems( elementSubRegion.fracturedElementsList().toViewConst()),
 
   61     m_cellsToEmbeddedSurfaces( elementSubRegion.embeddedSurfacesList().toViewConst() ),
 
   62     m_normalVector( embeddedSurfSubRegion.getNormalVector().toViewConst()),
 
   63     m_fracCenter( embeddedSurfSubRegion.getElementCenter().toViewConst()),
 
   64     m_fractureSurfaceArea( embeddedSurfSubRegion.getElementArea().toViewConst() ),
 
   65     m_connectivityIndex( embeddedSurfSubRegion.getConnectivityIndex() )
 
  103   template< 
typename POLICY,
 
  104             typename KERNEL_TYPE >
 
  109     forAll< POLICY >( kernelComponent.m_fracturedElems.size(),
 
  113       localIndex k = kernelComponent.m_fracturedElems[i];
 
  116       kernelComponent.setup( k, stack );
 
  118       real64 averageDistance = 0.0;
 
  125       kernelComponent.setConnectivityIndex( k, averageDistance );
 
  141       localIndex const localNodeIndex = m_elemsToNodes( k, a );
 
  143       for( 
int i=0; i<3; ++i )
 
  145         stack.
xLocal[ a ][ i ] = m_X[ localNodeIndex ][ i ];
 
  159                           real64 const (&point)[3] )
 const 
  161     localIndex const embSurfIndex = m_cellsToEmbeddedSurfaces[k][0];
 
  162     real64 pointToFracCenter[3];
 
  163     LvArray::tensorOps::copy< 3 >( pointToFracCenter, point );
 
  164     LvArray::tensorOps::subtract< 3 >( pointToFracCenter, m_fracCenter[embSurfIndex] );
 
  165     return LvArray::math::abs( LvArray::tensorOps::AiBi< 3 >( pointToFracCenter, m_normalVector[embSurfIndex] ));
 
  179     real64 parentSamplingPointCoord[3] = {0.0, 0.0, 0.0};
 
  180     FE_TYPE::getSamplingPointCoordInParentSpace( np, parentSamplingPointCoord );
 
  184     FE_TYPE::calcN( parentSamplingPointCoord, N );
 
  191       for( 
int i =0; i < 3; i++ )
 
  206                              real64 const averageDistance )
 const 
  208     localIndex const embSurfIndex = m_cellsToEmbeddedSurfaces[k][0];
 
  209     m_connectivityIndex[embSurfIndex] = 2. * m_fractureSurfaceArea[embSurfIndex] / averageDistance;
 
  215   FE_TYPE 
const & m_finiteElementSpace;
 
  224   traits::ViewTypeConst< typename CellElementSubRegion::NodeMapType::base_type > 
const m_elemsToNodes;
 
#define GEOS_HOST_DEVICE
Marks a host-device function.
 
#define GEOS_MARK_FUNCTION
Mark function with both Caliper and NVTX if enabled.
 
Kernel to compute EDFM connectivity index.
 
static constexpr int numNodesPerElem
number of nodes per element
 
GEOS_HOST_DEVICE void setConnectivityIndex(localIndex const k, real64 const averageDistance) const
Set the Connectivity Index.
 
static void launchCIComputationKernel(KERNEL_TYPE &kernelComponent)
launch of CI calculation
 
CIcomputationKernel(FE_TYPE const &finiteElementSpace, NodeManager const &nodeManager, CellElementSubRegion const &elementSubRegion, EmbeddedSurfaceSubRegion &embeddedSurfSubRegion)
Construct a new CIcomputationKernel object.
 
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
set up the kernel object by copying global values in the stack.
 
GEOS_HOST_DEVICE real64 computeDistance(localIndex const k, real64 const (&point)[3]) const
computes the distance between a point inside the element and the cut surface k.
 
static constexpr int numSamplingPoints
number of sampling points per element
 
GEOS_HOST_DEVICE void samplingPointCoord(integer const np, StackVariables &stack) const
computes coordinates of the sampling point in the physical space.
 
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
 
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).
 
ElementType
Denotes type of cell/element shape.
 
LvArray::SortedArrayView< T, localIndex, LvArray::ChaiBuffer > SortedArrayView
A sorted array view of local indices.
 
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
 
int integer
Signed integer type.
 
real64 xLocal[numNodesPerElem][3]
C-array stack storage for element local the nodal positions.
 
GEOS_HOST_DEVICE StackVariables()
Constructor.
 
real64 samplingPointCoord[3]
C-array stack storage for sampling points coordinates.