GEOS
CompositionalMultiphaseBase.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 Total, S.A
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
20 #ifndef GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONALMULTIPHASEBASE_HPP_
21 #define GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONALMULTIPHASEBASE_HPP_
22 
25 
26 namespace geos
27 {
28 
29 //START_SPHINX_INCLUDE_00
36 {
37 public:
38 
44  CompositionalMultiphaseBase( const string & name,
45  Group * const parent );
46 
49 
52 
55 
58 
61 
65  virtual ~CompositionalMultiphaseBase() override = default;
66 
67 //START_SPHINX_INCLUDE_01
68 
69  virtual void registerDataOnMesh( Group & meshBodies ) override;
70 
78  virtual void
79  implicitStepSetup( real64 const & time_n,
80  real64 const & dt,
81  DomainPartition & domain ) override;
82 
83  virtual void
84  assembleSystem( real64 const time_n,
85  real64 const dt,
86  DomainPartition & domain,
87  DofManager const & dofManager,
89  arrayView1d< real64 > const & localRhs ) override;
90 
91  virtual void
93  real64 const dt,
94  DomainPartition & domain,
95  DofManager const & dofManager,
97  arrayView1d< real64 > const & localRhs ) override;
98 
99  virtual void
101 
102  virtual void
104  real64 const & dt,
105  DomainPartition & domain ) override;
106 
112 
118 
123  void updateFluidModel( ObjectManagerBase & dataGroup ) const;
124 
129  void updateRelPermModel( ObjectManagerBase & dataGroup ) const;
130 
135  void updateCapPressureModel( ObjectManagerBase & dataGroup ) const;
136 
141  void updateCompAmount( ElementSubRegionBase & subRegion ) const;
142 
147  void updateEnergy( ElementSubRegionBase & subRegion ) const;
148 
154 
159  virtual void updatePhaseMobility( ObjectManagerBase & dataGroup ) const = 0;
160 
161  real64 updateFluidState( ElementSubRegionBase & subRegion ) const;
162 
163  virtual void saveConvergedState( ElementSubRegionBase & subRegion ) const override final;
164 
165  virtual void saveSequentialIterationState( DomainPartition & domain ) override final;
166 
167  virtual void updateState( DomainPartition & domain ) override final;
168 
174 
179  integer numFluidPhases() const { return m_numPhases; }
180 
186 
190  virtual units::Unit getMassUnit() const override
191  { return m_useMass ? units::Unit::Mass : units::Unit::Mole; }
192 
203  DofManager const & dofManager,
204  CRSMatrixView< real64, globalIndex const > const & localMatrix,
205  arrayView1d< real64 > const & localRhs ) const;
206 
216  virtual void
218  DomainPartition const & domain,
219  DofManager const & dofManager,
220  CRSMatrixView< real64, globalIndex const > const & localMatrix,
221  arrayView1d< real64 > const & localRhs ) const = 0;
222 
231  virtual void
233  DomainPartition const & domain,
234  DofManager const & dofManager,
235  CRSMatrixView< real64, globalIndex const > const & localMatrix,
236  arrayView1d< real64 > const & localRhs ) const = 0;
240  {
241  static constexpr char const * elemDofFieldString() { return "compositionalVariables"; }
242 
243  // inputs
244 
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"; }
250 
251 
252  // time stepping controls
253 
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"; }
260 
261 
262  // nonlinear solver parameters
263 
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"; }
274 
275  };
276 
285  void initializeFluidState( MeshLevel & mesh, DomainPartition & domain, arrayView1d< string const > const & regionNames );
286 
291 
301  void applyDirichletBC( real64 const time,
302  real64 const dt,
303  DofManager const & dofManager,
304  DomainPartition & domain,
305  CRSMatrixView< real64, globalIndex const > const & localMatrix,
306  arrayView1d< real64 > const & localRhs ) const;
307 
317  void applySourceFluxBC( real64 const time,
318  real64 const dt,
319  DofManager const & dofManager,
320  DomainPartition & domain,
321  CRSMatrixView< real64, globalIndex const > const & localMatrix,
322  arrayView1d< real64 > const & localRhs ) const;
323 
333  virtual void applyAquiferBC( real64 const time,
334  real64 const dt,
335  DofManager const & dofManager,
336  DomainPartition & domain,
337  CRSMatrixView< real64, globalIndex const > const & localMatrix,
338  arrayView1d< real64 > const & localRhs ) const = 0;
339 
352  real64 const dt,
353  DofManager const & dofManager,
354  DomainPartition & domain,
355  CRSMatrixView< real64, globalIndex const > const & localMatrix,
356  arrayView1d< real64 > const & localRhs ) const;
357 
358 
364 
365  virtual real64 setNextDtBasedOnStateChange( real64 const & currentDt,
366  DomainPartition & domain ) override;
367 
368  void computeCFLNumbers( DomainPartition & domain, real64 const & dt, real64 & maxPhaseCFL, real64 & maxCompCFL );
369 
376  real64 setNextDt( real64 const & currentDt,
377  DomainPartition & domain ) override;
378 
379  virtual real64 setNextDtBasedOnCFL( real64 const & currentDt,
380  DomainPartition & domain ) override;
381 
383 
384  integer useSimpleAccumulation() const { return m_useSimpleAccumulation; }
385 
386  integer useTotalMassEquation() const { return m_useTotalMassEquation; }
387 
388  virtual bool checkSequentialSolutionIncrements( DomainPartition & domain ) const override;
389 
390 protected:
391 
392  virtual void postInputInitialization() override;
393 
394  virtual void initializePreSubGroups() override;
395 
396 
403  void validateConstitutiveModels( DomainPartition const & domain ) const;
404 
409  void initializeAquiferBC( constitutive::ConstitutiveManager const & cm ) const;
410 
420  template< typename OBJECT_TYPE >
421  void applyFieldValue( real64 const & time_n,
422  real64 const & dt,
423  MeshLevel & mesh,
424  char const logMessage[],
425  string const fieldKey,
426  string const boundaryFieldKey ) const;
427 
430 
433 
436 
439 
442 
445 
448 
451 
454 
457 
460 
463 
466 
469 
472 
475 
478 
481 
484 
487 
490 
493  real64 m_maxSequentialCompDensChange;
494 
497 
498 private:
499 
505  bool validateDirichletBC( DomainPartition & domain,
506  real64 const time ) const;
507 
508  virtual void setConstitutiveNames( ElementSubRegionBase & subRegion ) const override;
509 
510 };
511 
512 template< typename OBJECT_TYPE >
514  real64 const & dt,
515  MeshLevel & mesh,
516  char const logMessage[],
517  string const fieldKey,
518  string const boundaryFieldKey ) const
519 {
521 
522  fsManager.apply< OBJECT_TYPE >( time_n + dt,
523  mesh,
524  fieldKey,
525  [&]( FieldSpecificationBase const & fs,
526  string const & setName,
528  OBJECT_TYPE & targetGroup,
529  string const & )
530  {
531  if( fs.getLogLevel() >= 1 && m_nonlinearSolverParameters.m_numNewtonIterations == 0 )
532  {
533  globalIndex const numTargetElems = MpiWrapper::sum< globalIndex >( lset.size() );
534  GEOS_LOG_RANK_0( GEOS_FMT( logMessage,
535  getName(), time_n+dt, fs.getCatalogName(), fs.getName(),
536  setName, targetGroup.getName(), fs.getScale(), numTargetElems ) );
537  }
538 
539  // Specify the bc value of the field
540  fs.applyFieldValue< FieldSpecificationEqual,
541  parallelDevicePolicy<> >( lset,
542  time_n + dt,
543  targetGroup,
544  boundaryFieldKey );
545  } );
546 }
547 
548 
549 } // namespace geos
550 
551 #endif //GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONALMULTIPHASEBASE_HPP_
#define GEOS_LOG_RANK_0(msg)
Log a message on screen on rank 0.
Definition: Logger.hpp:101
Unit
Enumerator of available unit types. Units are in SI by default.
Definition: Units.hpp:54
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 &currentDt, 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 &currentDt, 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 &currentDt, DomainPartition &domain) override
function to set the next dt based on state change
void initializeFluidState(MeshLevel &mesh, DomainPartition &domain, arrayView1d< string const > const &regionNames)
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,...
Definition: DofManager.hpp:44
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.
Definition: MeshLevel.hpp:42
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.
Definition: Group.hpp:1329
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.
Definition: DataTypes.hpp:180
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:310
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:88
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
LvArray::SortedArrayView< T, localIndex, LvArray::ChaiBuffer > SortedArrayView
A sorted array view of local indices.
Definition: DataTypes.hpp:271