20 #ifndef GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONALMULTIPHASEBASE_HPP_
21 #define GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONALMULTIPHASEBASE_HPP_
45 Group *
const parent );
191 {
return m_useMass ? units::Unit::Mass : units::Unit::Mole; }
241 static constexpr
char const * elemDofFieldString() {
return "compositionalVariables"; }
245 static constexpr
char const * useMassFlagString() {
return "useMass"; }
246 static constexpr
char const * relPermNamesString() {
return "relPermNames"; }
247 static constexpr
char const * capPressureNamesString() {
return "capPressureNames"; }
248 static constexpr
char const * diffusionNamesString() {
return "diffusionNames"; }
249 static constexpr
char const * dispersionNamesString() {
return "dispersionNames"; }
254 static constexpr
char const * solutionChangeScalingFactorString() {
return "solutionChangeScalingFactor"; }
255 static constexpr
char const * targetRelativePresChangeString() {
return "targetRelativePressureChangeInTimeStep"; }
256 static constexpr
char const * targetRelativeTempChangeString() {
return "targetRelativeTemperatureChangeInTimeStep"; }
257 static constexpr
char const * targetPhaseVolFracChangeString() {
return "targetPhaseVolFractionChangeInTimeStep"; }
258 static constexpr
char const * targetRelativeCompDensChangeString() {
return "targetRelativeCompDensChangeInTimeStep"; }
259 static constexpr
char const * targetFlowCFLString() {
return "targetFlowCFL"; }
264 static constexpr
char const * maxCompFracChangeString() {
return "maxCompFractionChange"; }
265 static constexpr
char const * maxRelativePresChangeString() {
return "maxRelativePressureChange"; }
266 static constexpr
char const * maxRelativeTempChangeString() {
return "maxRelativeTemperatureChange"; }
267 static constexpr
char const * maxRelativeCompDensChangeString() {
return "maxRelativeCompDensChange"; }
268 static constexpr
char const * allowLocalCompDensChoppingString() {
return "allowLocalCompDensityChopping"; }
269 static constexpr
char const * useTotalMassEquationString() {
return "useTotalMassEquation"; }
270 static constexpr
char const * useSimpleAccumulationString() {
return "useSimpleAccumulation"; }
271 static constexpr
char const * minCompDensString() {
return "minCompDens"; }
272 static constexpr
char const * maxSequentialCompDensChangeString() {
return "maxSequentialCompDensChange"; }
273 static constexpr
char const * minScalingFactorString() {
return "minScalingFactor"; }
420 template<
typename OBJECT_TYPE >
424 char const logMessage[],
425 string const fieldKey,
426 string const boundaryFieldKey )
const;
493 real64 m_maxSequentialCompDensChange;
506 real64 const time )
const;
512 template<
typename OBJECT_TYPE >
516 char const logMessage[],
517 string const fieldKey,
518 string const boundaryFieldKey )
const
522 fsManager.
apply< OBJECT_TYPE >( time_n + dt,
526 string const & setName,
528 OBJECT_TYPE & targetGroup,
533 globalIndex const numTargetElems = MpiWrapper::sum< globalIndex >( lset.size() );
535 getName(), time_n+dt, fs.getCatalogName(), fs.getName(),
536 setName, targetGroup.getName(), fs.getScale(), numTargetElems ) );
541 parallelDevicePolicy<> >( lset,
#define GEOS_LOG_RANK_0(msg)
Log a message on screen on rank 0.
Unit
Enumerator of available unit types. Units are in SI by default.
integer m_hasDispersion
flag to determine whether or not to apply dispersion
integer m_hasCapPressure
flag to determine whether or not to apply capillary pressure
void initializeAquiferBC(constitutive::ConstitutiveManager const &cm) const
Initialize the aquifer boundary condition (gravity vector, water phase index)
real64 m_solutionChangeScalingFactor
damping factor for solution change targets
CompositionalMultiphaseBase(CompositionalMultiphaseBase &&)=default
default move constructor
void applyFieldValue(real64 const &time_n, real64 const &dt, MeshLevel &mesh, char const logMessage[], string const fieldKey, string const boundaryFieldKey) const
Utility function that encapsulates the call to FieldSpecificationBase::applyFieldValue in BC applicat...
virtual void initializePreSubGroups() override
Called by Initialize() prior to initializing sub-Groups.
real64 m_targetPhaseVolFracChange
target (absolute) change in phase volume fraction in a time step
virtual real64 setNextDtBasedOnCFL(real64 const ¤tDt, DomainPartition &domain) override
function to set the next dt based on state change
virtual void postInputInitialization() override
void keepVariablesConstantDuringInitStep(real64 const time, real64 const dt, DofManager const &dofManager, DomainPartition &domain, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const
Function to fix the initial state during the initialization step in coupled problems.
integer m_useMass
flag indicating whether mass or molar formulation should be used
integer m_useSimpleAccumulation
flag indicating whether simple accumulation form is used
real64 m_maxRelativeTempChange
maximum (relative) change in temperature in a Newton iteration
virtual bool checkSequentialSolutionIncrements(DomainPartition &domain) const override
Check if the solution increments are ok to use.
real64 m_maxRelativePresChange
maximum (relative) change in pressure in a Newton iteration
real64 m_targetRelativePresChange
target (relative) change in pressure in a time step
void validateConstitutiveModels(DomainPartition const &domain) const
Utility function that checks the consistency of the constitutive models.
integer m_hasDiffusion
flag to determine whether or not to apply diffusion
real64 m_minScalingFactor
minimum value of the scaling factor obtained by enforcing maxCompFracChange
virtual void applyAquiferBC(real64 const time, real64 const dt, DofManager const &dofManager, DomainPartition &domain, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const =0
Apply aquifer boundary conditions to the system.
virtual void registerDataOnMesh(Group &meshBodies) override
Register wrappers that contain data on the mesh objects.
void applySourceFluxBC(real64 const time, real64 const dt, DofManager const &dofManager, DomainPartition &domain, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const
Apply source flux boundary conditions to the system.
real64 setNextDt(real64 const ¤tDt, DomainPartition &domain) override
function to set the next time step size
virtual void initializePostInitialConditionsPreSubGroups() override
Called by InitializePostInitialConditions() prior to initializing sub-Groups.
real64 m_maxCompFracChange
maximum (absolute) change in a component fraction in a Newton iteration
real64 m_sequentialCompDensChange
maximum (absolute) component density change in a sequential iteration
real64 m_targetRelativeCompDensChange
target (relative) change in component density in a time step
virtual real64 setNextDtBasedOnStateChange(real64 const ¤tDt, DomainPartition &domain) override
function to set the next dt based on state change
void initializeFluidState(MeshLevel &mesh, DomainPartition &domain, arrayView1d< string const > const ®ionNames)
Initialize all variables from initial conditions.
CompositionalMultiphaseBase()=delete
deleted default constructor
integer m_allowCompDensChopping
flag indicating whether local (cell-wise) chopping of negative compositions is allowed
string m_referenceFluidModelName
name of the fluid constitutive model used as a reference for component/phase description
integer m_numComponents
the number of fluid components
CompositionalMultiphaseBase(CompositionalMultiphaseBase const &)=delete
deleted copy constructor
integer m_useTotalMassEquation
flag indicating whether total mass equation is used
virtual ~CompositionalMultiphaseBase() override=default
default destructor
real64 m_targetFlowCFL
the targeted CFL for timestep
void chopNegativeDensities(DomainPartition &domain)
Sets all the negative component densities (if any) to zero.
real64 m_maxRelativeCompDensChange
maximum (relative) change in component density in a Newton iteration
real64 m_minCompDens
minimum allowed global component density
CompositionalMultiphaseBase(const string &name, Group *const parent)
main constructor for Group Objects
integer m_numPhases
the max number of fluid phases
real64 m_targetRelativeTempChange
target (relative) change in temperature in a time step
void applyDirichletBC(real64 const time, real64 const dt, DofManager const &dofManager, DomainPartition &domain, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const
Function to perform the Application of Dirichlet type BC's.
CompositionalMultiphaseBase & operator=(CompositionalMultiphaseBase const &)=delete
deleted assignment operator
void computeHydrostaticEquilibrium()
Compute the hydrostatic equilibrium using the compositions and temperature input tables.
CompositionalMultiphaseBase & operator=(CompositionalMultiphaseBase &&)=delete
deleted move operator
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
void apply(real64 const time, MeshLevel &mesh, string const &fieldName, LAMBDA &&lambda) const
This function is the main driver for the field applications.
static FieldSpecificationManager & getInstance()
Class facilitating the representation of a multi-level discretization of a MeshBody.
integer m_numNewtonIterations
The number of nonlinear iterations that have been exectued.
The ObjectManagerBase is the base object of all object managers in the mesh data hierachy.
NonlinearSolverParameters m_nonlinearSolverParameters
Nonlinear solver parameters.
string const & getName() const
Get group name.
virtual void resetStateToBeginningOfStep(DomainPartition &domain) override
reset state of physics back to the beginning of the step.
void assembleAccumulationAndVolumeBalanceTerms(DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const
assembles the accumulation and volume balance terms for all cells
virtual void implicitStepSetup(real64 const &time_n, real64 const &dt, DomainPartition &domain) override
function to perform setup for implicit timestep
real64 updatePhaseVolumeFraction(ObjectManagerBase &dataGroup) const
Recompute phase volume fractions (saturations) from constitutive and primary variables.
virtual units::Unit getMassUnit() const override
virtual void assembleSystem(real64 const time_n, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) override
function to assemble the linear system matrix and rhs
virtual void assembleStabilizedFluxTerms(real64 const dt, DomainPartition const &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const =0
assembles the flux terms for all cells with pressure jump stabilization
virtual void saveSequentialIterationState(DomainPartition &domain) override final
Utility function to save the iteration state (useful for sequential simulations)
virtual void implicitStepComplete(real64 const &time, real64 const &dt, DomainPartition &domain) override
perform cleanup for implicit timestep
integer numFluidComponents() const
Getter for the number of fluid components (species)
virtual void applyBoundaryConditions(real64 const time_n, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) override
apply boundary condition to system
virtual void assembleFluxTerms(real64 const dt, DomainPartition const &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const =0
assembles the flux terms for all cells
void updateFluidModel(ObjectManagerBase &dataGroup) const
Update all relevant fluid models using current values of pressure and composition.
virtual void updatePhaseMobility(ObjectManagerBase &dataGroup) const =0
Recompute phase mobility from constitutive and primary variables.
void updateRelPermModel(ObjectManagerBase &dataGroup) const
Update all relevant relperm models using current values of phase volume fraction.
virtual void updateState(DomainPartition &domain) override final
Recompute all dependent quantities from primary variables (including constitutive models)
void updateCapPressureModel(ObjectManagerBase &dataGroup) const
Update all relevant capillary pressure models using current values of phase volume fraction.
void updateEnergy(ElementSubRegionBase &subRegion) const
Update energy.
void updateSolidInternalEnergyModel(ObjectManagerBase &dataGroup) const
Update all relevant solid internal energy models using current values of temperature.
string referenceFluidModelName() const
Getter for the name of the reference fluid model name.
virtual void saveConvergedState(ElementSubRegionBase &subRegion) const override final
Utility function to save the converged state.
void updateCompAmount(ElementSubRegionBase &subRegion) const
Update components mass/moles.
integer numFluidPhases() const
Getter for the number of fluid phases.
void updateGlobalComponentFraction(ObjectManagerBase &dataGroup) const
Recompute global component fractions from primary variables (component densities)
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.
std::int32_t integer
Signed integer type.
LvArray::SortedArrayView< T, localIndex, LvArray::ChaiBuffer > SortedArrayView
A sorted array view of local indices.