20 #ifndef GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONAL_FLUXCOMPUTEZFORMULATIONKERNEL_HPP
21 #define GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONAL_FLUXCOMPUTEZFORMULATIONKERNEL_HPP
25 #include "codingUtilities/Utilities.hpp"
28 #include "common/GEOS_RAJA_Interface.hpp"
29 #include "constitutive/fluid/multifluid/MultiFluidFields.hpp"
34 #include "physicsSolvers/fluidFlow/kernels/compositional/KernelLaunchSelectors.hpp"
43 namespace isothermalCompositionalMultiphaseFVMKernels
53 template<
integer NUM_COMP,
integer NUM_DOF,
typename STENCILWRAPPER >
96 STENCILWRAPPER
const & stencilWrapper,
97 DofNumberAccessor
const & dofNumberAccessor,
105 BitFlags< KernelFlags > kernelFlags )
115 m_permeability( permeabilityAccessors.get( fields::permeability::permeability {} ) ),
116 m_dPerm_dPres( permeabilityAccessors.get( fields::permeability::dPerm_dPressure {} ) ),
117 m_phaseMob( compFlowAccessors.get( fields::flow::phaseMobility {} ) ),
118 m_dPhaseMob( compFlowAccessors.get( fields::flow::dPhaseMobility {} ) ),
119 m_phaseMassDens( multiFluidAccessors.get( fields::multifluid::phaseMassDensity {} ) ),
120 m_dPhaseMassDens( multiFluidAccessors.get( fields::multifluid::dPhaseMassDensity {} ) ),
121 m_phaseCapPressure( capPressureAccessors.get( fields::cappres::phaseCapPressure {} ) ),
122 m_dPhaseCapPressure_dPhaseVolFrac( capPressureAccessors.get( fields::cappres::dPhaseCapPressure_dPhaseVolFraction {} ) ),
124 m_seri( stencilWrapper.getElementRegionIndices() ),
125 m_sesri( stencilWrapper.getElementSubRegionIndices() ),
126 m_sei( stencilWrapper.getElementIndices() )
225 template<
typename FUNC = NoOpFunc >
230 FUNC && compFluxKernelOp = NoOpFunc{} )
const
275 isothermalCompositionalMultiphaseFVMKernelUtilities::PPUPhaseFluxZFormulation::compute< numComp, numFluxSupportPoints >
278 m_kernelFlags.isSet( KernelFlags::CapPressure ),
279 m_kernelFlags.isSet( KernelFlags::CheckPhasePresenceInGravity ),
301 compFluxKernelOp( ip, k, seri, sesri, sei, connectionIndex,
302 k_up, seri[k_up], sesri[k_up], sei[k_up], potGrad,
303 phaseFlux, dPhaseFlux_dP, dPhaseFlux_dC );
324 localIndex const localDofIndexComp = localDofIndexPres + jc + 1;
341 template<
typename FUNC = NoOpFunc >
346 FUNC && assemblyKernelOp = NoOpFunc{} )
const
348 using namespace compositionalMultiphaseUtilities;
350 if( m_kernelFlags.isSet( KernelFlags::TotalMassEquation ) )
365 if(
m_ghostRank[
m_seri( iconn, i )][m_sesri( iconn, i )][m_sei( iconn, i )] < 0 )
374 RAJA::atomicAdd( parallelDeviceAtomic{}, &
m_localRhs[localRow + ic],
376 m_localMatrix.addToRowBinarySearchUnsorted< parallelDeviceAtomic >
384 assemblyKernelOp( i, localRow );
396 template<
typename POLICY,
typename KERNEL_TYPE >
399 KERNEL_TYPE
const & kernelComponent )
404 typename KERNEL_TYPE::StackVariables stack( kernelComponent.stencilSize( iconn ),
405 kernelComponent.numPointsInFlux( iconn ) );
407 kernelComponent.setup( iconn, stack );
408 kernelComponent.computeFlux( iconn, stack );
409 kernelComponent.complete( iconn, stack );
437 typename STENCILWRAPPER::IndexContainerViewConstType
const m_seri;
438 typename STENCILWRAPPER::IndexContainerViewConstType
const m_sesri;
439 typename STENCILWRAPPER::IndexContainerViewConstType
const m_sei;
466 template<
typename POLICY,
typename STENCILWRAPPER >
471 string const & dofKey,
472 BitFlags< KernelFlags > kernelFlags,
473 string const & solverName,
475 STENCILWRAPPER
const & stencilWrapper,
480 isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComps, [&](
auto NC )
482 integer constexpr NUM_COMP = NC();
483 integer constexpr NUM_DOF = NC() + 1;
487 dofNumberAccessor.setName( solverName +
"/accessors/" + dofKey );
490 typename kernelType::CompFlowAccessors compFlowAccessors( elemManager, solverName );
491 typename kernelType::MultiFluidAccessors multiFluidAccessors( elemManager, solverName );
492 typename kernelType::CapPressureAccessors capPressureAccessors( elemManager, solverName );
493 typename kernelType::PermeabilityAccessors permeabilityAccessors( elemManager, solverName );
495 kernelType kernel( numPhases, rankOffset, stencilWrapper, dofNumberAccessor,
496 compFlowAccessors, multiFluidAccessors, capPressureAccessors, permeabilityAccessors,
497 dt, localMatrix, localRhs, kernelFlags );
498 kernelType::template launch< POLICY >( stencilWrapper.size(), kernel );
#define GEOS_HOST_DEVICE
Marks a host-device function.
#define GEOS_ASSERT_GT(lhs, rhs)
Assert that one value compares greater than the other in debug builds.
#define GEOS_ASSERT_GE(lhs, rhs)
Assert that one value compares greater than or equal to the other in debug builds.
#define GEOS_MARK_FUNCTION
Mark function with both Caliper and NVTX if enabled.
The ElementRegionManager class provides an interface to ObjectManagerBase in order to manage ElementR...
ElementViewAccessor< ArrayView< T const, NDIM, getUSD< PERM > > > constructArrayViewAccessor(string const &name, string const &neighborName=string()) const
This is a function to construct a ElementViewAccessor to access array data registered on the mesh.
array1d< array1d< VIEWTYPE > > ElementViewAccessor
The ElementViewAccessor at the ElementRegionManager level is an array of array of VIEWTYPE.
A struct to automatically construct and store element view accessors.
A struct to automatically construct and store element view accessors.
Base class for FluxComputeKernel that holds all data not dependent on template parameters (like stenc...
ElementViewConst< arrayView2d< real64 const, compflow::USD_PHASE > > const m_phaseVolFrac
Views on phase volume fractions.
real64 const m_dt
Time step size.
ElementViewConst< arrayView1d< globalIndex const > > const m_dofNumber
Views on dof numbers.
ElementViewConst< arrayView1d< real64 const > > const m_pres
Views on pressure.
ElementViewConst< arrayView1d< integer const > > const m_ghostRank
Views on ghost rank numbers and gravity coefficients.
arrayView1d< real64 > const m_localRhs
View on the local RHS.
ElementViewConst< arrayView4d< real64 const, constitutive::multifluid::USD_PHASE_COMP > > const m_phaseCompFrac
Views on phase component fractions.
globalIndex const m_rankOffset
Offset for my MPI rank.
CRSMatrixView< real64, globalIndex const > const m_localMatrix
View on the local CRS matrix.
ElementRegionManager::ElementViewConst< VIEWTYPE > ElementViewConst
The type for element-based data. Consists entirely of ArrayView's.
integer const m_numPhases
Number of fluid phases.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
StackArray< T, 2, MAXSIZE > stackArray2d
Alias for 2D stack array.
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
StackArray< T, 1, MAXSIZE > stackArray1d
Alias for 1D stack array.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
std::int32_t integer
Signed integer type.