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 TotalEnergies
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 
71  virtual void registerDataForCFL( Group & meshBodies ) { GEOS_UNUSED_VAR( meshBodies ); }
72 
80  virtual void
81  implicitStepSetup( real64 const & time_n,
82  real64 const & dt,
83  DomainPartition & domain ) override;
84 
85  virtual void
86  assembleSystem( real64 const time_n,
87  real64 const dt,
88  DomainPartition & domain,
89  DofManager const & dofManager,
91  arrayView1d< real64 > const & localRhs ) override;
92 
93  virtual void
95  real64 const dt,
96  DomainPartition & domain,
97  DofManager const & dofManager,
99  arrayView1d< real64 > const & localRhs ) override;
100 
101  virtual void
103 
104  virtual void
106  real64 const & dt,
107  DomainPartition & domain ) override;
108 
114 
120 
125  void updateFluidModel( ObjectManagerBase & dataGroup ) const;
126 
131  void updateRelPermModel( ObjectManagerBase & dataGroup ) const;
132 
137  void updateCapPressureModel( ObjectManagerBase & dataGroup ) const;
138 
143  void updateCompAmount( ElementSubRegionBase & subRegion ) const;
144 
149  void updateEnergy( ElementSubRegionBase & subRegion ) const;
150 
156 
161  virtual void updatePhaseMobility( ObjectManagerBase & dataGroup ) const = 0;
162 
163  real64 updateFluidState( ElementSubRegionBase & subRegion ) const;
164 
165  virtual void saveConvergedState( ElementSubRegionBase & subRegion ) const override final;
166 
167  virtual void saveSequentialIterationState( DomainPartition & domain ) override final;
168 
169  virtual void updateState( DomainPartition & domain ) override final;
170 
176 
181  integer numFluidPhases() const { return m_numPhases; }
182 
188 
192  virtual units::Unit getMassUnit() const override
193  { return m_useMass ? units::Unit::Mass : units::Unit::Mole; }
194 
205  DofManager const & dofManager,
206  CRSMatrixView< real64, globalIndex const > const & localMatrix,
207  arrayView1d< real64 > const & localRhs ) const;
208 
218  virtual void
220  DomainPartition const & domain,
221  DofManager const & dofManager,
222  CRSMatrixView< real64, globalIndex const > const & localMatrix,
223  arrayView1d< real64 > const & localRhs ) const = 0;
224 
233  virtual void
235  DomainPartition const & domain,
236  DofManager const & dofManager,
237  CRSMatrixView< real64, globalIndex const > const & localMatrix,
238  arrayView1d< real64 > const & localRhs ) const = 0;
242  {
243  static constexpr char const * elemDofFieldString() { return "compositionalVariables"; }
244 
245  // inputs
246 
247  static constexpr char const * useMassFlagString() { return "useMass"; }
248  static constexpr char const * relPermNamesString() { return "relPermNames"; }
249  static constexpr char const * capPressureNamesString() { return "capPressureNames"; }
250  static constexpr char const * diffusionNamesString() { return "diffusionNames"; }
251  static constexpr char const * dispersionNamesString() { return "dispersionNames"; }
252 
253 
254  // time stepping controls
255 
256  static constexpr char const * solutionChangeScalingFactorString() { return "solutionChangeScalingFactor"; }
257  static constexpr char const * targetRelativePresChangeString() { return "targetRelativePressureChangeInTimeStep"; }
258  static constexpr char const * targetRelativeTempChangeString() { return "targetRelativeTemperatureChangeInTimeStep"; }
259  static constexpr char const * targetPhaseVolFracChangeString() { return "targetPhaseVolFractionChangeInTimeStep"; }
260  static constexpr char const * targetRelativeCompDensChangeString() { return "targetRelativeCompDensChangeInTimeStep"; }
261  static constexpr char const * targetFlowCFLString() { return "targetFlowCFL"; }
262 
263 
264  // nonlinear solver parameters
265 
266  static constexpr char const * maxCompFracChangeString() { return "maxCompFractionChange"; }
267  static constexpr char const * maxRelativePresChangeString() { return "maxRelativePressureChange"; }
268  static constexpr char const * maxRelativeTempChangeString() { return "maxRelativeTemperatureChange"; }
269  static constexpr char const * maxRelativeCompDensChangeString() { return "maxRelativeCompDensChange"; }
270  static constexpr char const * allowLocalCompDensChoppingString() { return "allowLocalCompDensityChopping"; }
271  static constexpr char const * useTotalMassEquationString() { return "useTotalMassEquation"; }
272  static constexpr char const * useSimpleAccumulationString() { return "useSimpleAccumulation"; }
273  static constexpr char const * minCompDensString() { return "minCompDens"; }
274  static constexpr char const * maxSequentialCompDensChangeString() { return "maxSequentialCompDensChange"; }
275  static constexpr char const * minScalingFactorString() { return "minScalingFactor"; }
276 
277  };
278 
287  virtual void initializeFluidState( MeshLevel & mesh, arrayView1d< string const > const & regionNames ) override;
288 
289  virtual void initializeThermalState( MeshLevel & mesh, arrayView1d< string const > const & regionNames ) override;
290 
294  virtual void computeHydrostaticEquilibrium( DomainPartition & domain ) override;
295 
305  void applyDirichletBC( real64 const time,
306  real64 const dt,
307  DofManager const & dofManager,
308  DomainPartition & domain,
309  CRSMatrixView< real64, globalIndex const > const & localMatrix,
310  arrayView1d< real64 > const & localRhs ) const;
311 
321  void applySourceFluxBC( real64 const time,
322  real64 const dt,
323  DofManager const & dofManager,
324  DomainPartition & domain,
325  CRSMatrixView< real64, globalIndex const > const & localMatrix,
326  arrayView1d< real64 > const & localRhs ) const;
327 
337  virtual void applyAquiferBC( real64 const time,
338  real64 const dt,
339  DofManager const & dofManager,
340  DomainPartition & domain,
341  CRSMatrixView< real64, globalIndex const > const & localMatrix,
342  arrayView1d< real64 > const & localRhs ) const = 0;
343 
356  real64 const dt,
357  DofManager const & dofManager,
358  DomainPartition & domain,
359  CRSMatrixView< real64, globalIndex const > const & localMatrix,
360  arrayView1d< real64 > const & localRhs ) const;
361 
362 
368 
369  void chopNegativeDensities( ElementSubRegionBase & subRegion );
370 
371  virtual real64 setNextDtBasedOnStateChange( real64 const & currentDt,
372  DomainPartition & domain ) override;
373 
374 
375  virtual void computeCFLNumbers( DomainPartition & domain, real64 const & dt, real64 & maxPhaseCFL, real64 & maxCompCFL )
376  {
377  GEOS_UNUSED_VAR( domain, dt, maxPhaseCFL, maxCompCFL );
378  GEOS_ERROR( GEOS_FMT( "{}: computeCFLNumbers is not implemented for {}", getDataContext(), getCatalogName() ) );
379  }
380 
382 
383  integer useSimpleAccumulation() const { return m_useSimpleAccumulation; }
384 
385  integer useTotalMassEquation() const { return m_useTotalMassEquation; }
386 
387  virtual bool checkSequentialSolutionIncrements( DomainPartition & domain ) const override;
388 
389 protected:
390 
391  virtual void postInputInitialization() override;
392 
393  virtual void initializePreSubGroups() override;
394 
395 
402  void validateConstitutiveModels( DomainPartition const & domain ) const;
403 
408  void initializeAquiferBC( constitutive::ConstitutiveManager const & cm ) const;
409 
419  template< typename OBJECT_TYPE >
420  void applyFieldValue( real64 const & time_n,
421  real64 const & dt,
422  MeshLevel & mesh,
423  char const logMessage[],
424  string const fieldKey,
425  string const boundaryFieldKey ) const;
426 
429 
432 
435 
438 
441 
444 
447 
450 
453 
456 
459 
462 
465 
468 
471 
474 
477 
480 
483 
486 
489 
492  real64 m_maxSequentialCompDensChange;
493 
494 private:
495 
501  bool validateDirichletBC( DomainPartition & domain,
502  real64 const time ) const;
503 
504  virtual void setConstitutiveNames( ElementSubRegionBase & subRegion ) const override;
505 
506 };
507 
508 template< typename OBJECT_TYPE >
510  real64 const & dt,
511  MeshLevel & mesh,
512  char const logMessage[],
513  string const fieldKey,
514  string const boundaryFieldKey ) const
515 {
517 
518  fsManager.apply< OBJECT_TYPE >( time_n + dt,
519  mesh,
520  fieldKey,
521  [&]( FieldSpecificationBase const & fs,
522  string const & setName,
524  OBJECT_TYPE & targetGroup,
525  string const & )
526  {
527  if( fs.getLogLevel() >= 1 && m_nonlinearSolverParameters.m_numNewtonIterations == 0 )
528  {
529  globalIndex const numTargetElems = MpiWrapper::sum< globalIndex >( lset.size() );
530  GEOS_LOG_RANK_0( GEOS_FMT( logMessage,
531  getName(), time_n+dt, fs.getCatalogName(), fs.getName(),
532  setName, targetGroup.getName(), fs.getScale(), numTargetElems ) );
533  }
534 
535  // Specify the bc value of the field
536  fs.applyFieldValue< FieldSpecificationEqual,
537  parallelDevicePolicy<> >( lset,
538  time_n + dt,
539  targetGroup,
540  boundaryFieldKey );
541  } );
542 }
543 
544 
545 } // namespace geos
546 
547 #endif //GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONALMULTIPHASEBASE_HPP_
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
Definition: Logger.hpp:157
#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 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.
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
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
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
virtual void computeHydrostaticEquilibrium(DomainPartition &domain) override
Compute the hydrostatic equilibrium using the compositions and temperature input tables.
real64 m_minCompDens
minimum allowed global component density
virtual void initializeFluidState(MeshLevel &mesh, arrayView1d< string const > const &regionNames) override
Initialize all variables from initial conditions.
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
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.
virtual string getCatalogName() const =0
NonlinearSolverParameters m_nonlinearSolverParameters
Nonlinear solver parameters.
DataContext const & getDataContext() const
Definition: Group.hpp:1343
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