20 #ifndef GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONAL_STABILIZEDFLUXCOMPUTEKERNEL_HPP
21 #define GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONAL_STABILIZEDFLUXCOMPUTEKERNEL_HPP
24 #include "constitutive/relativePermeability/RelativePermeabilityBase.hpp"
25 #include "constitutive/relativePermeability/RelativePermeabilityFields.hpp"
30 namespace stabilizedCompositionalMultiphaseFVMKernels
42 template<
integer NUM_COMP,
integer NUM_DOF,
typename STENCILWRAPPER >
53 template<
typename VIEWTYPE >
57 using DofNumberAccessor = AbstractBase::DofNumberAccessor;
65 fields::flow::elementStabConstant,
66 fields::flow::pressure_n >;
70 fields::multifluid::phaseDensity_n,
71 fields::multifluid::phaseCompFraction_n >;
78 using AbstractBase::m_kernelFlags;
82 using AbstractBase::m_gravCoef;
96 using Base::m_dPhaseMassDens;
122 STENCILWRAPPER
const & stencilWrapper,
123 DofNumberAccessor
const & dofNumberAccessor,
134 BitFlags< isothermalCompositionalMultiphaseFVMKernels::KernelFlags > kernelFlags )
141 capPressureAccessors,
142 permeabilityAccessors,
147 m_pres_n( stabCompFlowAccessors.get( fields::flow::pressure_n {} ) ),
148 m_phaseDens_n( stabMultiFluidAccessors.get( fields::multifluid::phaseDensity_n {} ) ),
149 m_phaseCompFrac_n( stabMultiFluidAccessors.get( fields::multifluid::phaseCompFraction_n {} ) ),
150 m_phaseRelPerm_n( relPermAccessors.get( fields::relperm::phaseRelPerm_n {} ) ),
152 m_elementStabConstant( stabCompFlowAccessors.get( fields::flow::elementStabConstant {} ) )
208 real64 const (&dPhaseFlux_dP)[2],
211 GEOS_UNUSED_VAR( k_up, potGrad, phaseFlux, dPhaseFlux_dP, dPhaseFlux_dC, er_up, esr_up, ei_up );
222 real64 dPresGradStab = 0.0;
236 m_elementStabConstant[er][esr][ei] * stabTrans[ke] * (
m_pres[er][esr][ei] -
m_pres_n[er][esr][ei]);
240 integer const k_up_stab = (dPresGradStab >= 0) ? 0 : 1;
242 localIndex const er_up_stab = seri[k_up_stab];
243 localIndex const esr_up_stab = sesri[k_up_stab];
246 bool const areInSameMacroElement = stencilMacroElements[0] == stencilMacroElements[1];
247 bool const isStabilizationActive = stencilMacroElements[0] >= 0 && stencilMacroElements[1] >= 0;
248 if( isStabilizationActive && areInSameMacroElement )
253 real64 const laggedUpwindCoef = m_phaseDens_n[er_up_stab][esr_up_stab][ei_up_stab][0][ip]
254 * m_phaseCompFrac_n[er_up_stab][esr_up_stab][ei_up_stab][0][ip][ic]
255 * m_phaseRelPerm_n[er_up_stab][esr_up_stab][ei_up_stab][0][ip];
256 stabFlux[ic] += dPresGradStab * laggedUpwindCoef;
260 real64 const tauStab = m_elementStabConstant[seri[ke]][sesri[ke]][sei[ke]];
261 dStabFlux_dP[ke][ic] += tauStab * stabTrans[ke] * laggedUpwindCoef;
272 stack.
localFlux[eqIndex0] += stabFlux[ic];
273 stack.
localFlux[eqIndex1] += -stabFlux[ic];
324 template<
typename POLICY,
typename STENCILWRAPPER >
329 string const & dofKey,
330 BitFlags< isothermalCompositionalMultiphaseFVMKernels::KernelFlags > kernelFlags,
331 string const & solverName,
333 STENCILWRAPPER
const & stencilWrapper,
338 isothermalCompositionalMultiphaseBaseKernels::
339 internal::kernelLaunchSelectorCompSwitch( numComps, [&](
auto NC )
341 integer constexpr NUM_COMP = NC();
342 integer constexpr NUM_DOF = NC() + 1;
346 dofNumberAccessor.setName( solverName +
"/accessors/" + dofKey );
349 typename KERNEL_TYPE::CompFlowAccessors compFlowAccessors( elemManager, solverName );
350 typename KERNEL_TYPE::MultiFluidAccessors multiFluidAccessors( elemManager, solverName );
351 typename KERNEL_TYPE::StabCompFlowAccessors stabCompFlowAccessors( elemManager, solverName );
352 typename KERNEL_TYPE::StabMultiFluidAccessors stabMultiFluidAccessors( elemManager, solverName );
353 typename KERNEL_TYPE::CapPressureAccessors capPressureAccessors( elemManager, solverName );
354 typename KERNEL_TYPE::PermeabilityAccessors permeabilityAccessors( elemManager, solverName );
355 typename KERNEL_TYPE::RelPermAccessors relPermAccessors( elemManager, solverName );
357 KERNEL_TYPE kernel( numPhases, rankOffset, stencilWrapper, dofNumberAccessor,
358 compFlowAccessors, stabCompFlowAccessors, multiFluidAccessors, stabMultiFluidAccessors,
359 capPressureAccessors, permeabilityAccessors, relPermAccessors,
360 dt, localMatrix, localRhs, kernelFlags );
361 KERNEL_TYPE::template launch< POLICY >( stencilWrapper.size(), kernel );
#define GEOS_HOST_DEVICE
Marks a host-device function.
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
#define GEOS_UNUSED_PARAM(X)
Mark an unused argument and silence compiler warnings.
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.
typename ElementViewAccessor< VIEWTYPE >::NestedViewTypeConst ElementViewConst
The ElementViewAccessor at the ElementRegionManager level is the type resulting from ElementViewAcces...
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...
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.
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.
ElementViewConst< arrayView3d< real64 const, compflow::USD_COMP_DC > > const m_dCompFrac_dCompDens
Views on derivatives of comp fractions.
integer const m_numPhases
Number of fluid phases.
Define the interface for the assembly kernel in charge of flux terms.
static constexpr localIndex maxStencilSize
Maximum number of points in the stencil.
GEOS_HOST_DEVICE void computeFlux(localIndex const iconn, StackVariables &stack, FUNC &&compFluxKernelOp=NoOpFunc{}) const
Compute the local flux contributions to the residual and Jacobian.
static constexpr integer numFluxSupportPoints
Number of flux support points (hard-coded for TFPA)
static constexpr integer numEqn
Compute time value for the number of equations (all of them, except the volume balance equation)
static constexpr integer numComp
Compile time value for the number of components.
ElementViewConst< arrayView2d< real64 const, compflow::USD_PHASE > > const m_phaseMob
Views on phase mobilities.
STENCILWRAPPER::IndexContainerViewConstType const m_seri
Connection to element maps.
static constexpr localIndex maxNumElems
Maximum number of elements at the face.
static constexpr integer numDof
Compute time value for the number of degrees of freedom.
STENCILWRAPPER const m_stencilWrapper
Reference to the stencil wrapper.
static constexpr localIndex maxNumConns
Maximum number of connections at the face.
static void createAndLaunch(integer const numComps, integer const numPhases, globalIndex const rankOffset, string const &dofKey, BitFlags< isothermalCompositionalMultiphaseFVMKernels::KernelFlags > kernelFlags, string const &solverName, ElementRegionManager const &elemManager, STENCILWRAPPER const &stencilWrapper, real64 const dt, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)
Create a new kernel and launch.
Define the interface for the assembly kernel in charge of flux terms.
ElementViewConst< arrayView1d< real64 const > > const m_pres_n
Views on flow properties at the previous converged time step.
ElementViewConst< arrayView1d< integer const > > const m_macroElementIndex
Views on the macroelement indices and stab constant.
ElementViewConst< arrayView1d< real64 const > > const m_pres
Views on pressure.
static constexpr integer numFluxSupportPoints
Number of flux support points (hard-coded for TFPA)
static constexpr integer numEqn
Compute time value for the number of equations (all of them, except the volume balance equation)
static constexpr integer numComp
Compile time value for the number of components.
GEOS_HOST_DEVICE void computeFlux(localIndex const iconn, StackVariables &stack) const
Compute the local flux contributions to the residual and Jacobian, including stabilization.
FluxComputeKernel(integer const numPhases, globalIndex const rankOffset, STENCILWRAPPER const &stencilWrapper, DofNumberAccessor const &dofNumberAccessor, CompFlowAccessors const &compFlowAccessors, StabCompFlowAccessors const &stabCompFlowAccessors, MultiFluidAccessors const &multiFluidAccessors, StabMultiFluidAccessors const &stabMultiFluidAccessors, CapPressureAccessors const &capPressureAccessors, PermeabilityAccessors const &permeabilityAccessors, RelPermAccessors const &relPermAccessors, real64 const &dt, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs, BitFlags< isothermalCompositionalMultiphaseFVMKernels::KernelFlags > kernelFlags)
Constructor for the kernel interface.
ElementRegionManager::ElementViewConst< VIEWTYPE > ElementViewConst
The type for element-based data. Consists entirely of ArrayView's.
static constexpr integer numDof
Compute time value for the number of degrees of freedom.
STENCILWRAPPER const m_stencilWrapper
Reference to the stencil wrapper.
static constexpr localIndex maxNumConns
Maximum number of connections at the face.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
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).
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.
Kernel variables (dof numbers, jacobian and residual) located on the stack.
real64 dTrans_dPres[maxNumConns][numFluxSupportPoints]
Derivatives of transmissibility with respect to pressure.
stackArray1d< globalIndex, maxNumElems *numDof > dofColIndices
Indices of the matrix rows/columns corresponding to the dofs in this face.
stackArray2d< real64, maxNumElems *numEqn *maxStencilSize *numDof > localFluxJacobian
Storage for the face local Jacobian matrix.
stackArray1d< real64, maxNumElems *numEqn > localFlux
Storage for the face local residual vector (all equations except volume balance)
localIndex const numConnectedElems
Number of elements connected at a given connection.
real64 transmissibility[maxNumConns][numFluxSupportPoints]
Transmissibility.
localIndex const stencilSize
Stencil size for a given connection.
real64 stabTransmissibility[maxNumConns][numFluxSupportPoints]
Stabilization transmissibility.