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 
24 #include "common/DataLayouts.hpp"
25 #include "constitutive/fluid/multifluid/Layouts.hpp"
26 #include "constitutive/relativePermeability/Layouts.hpp"
27 #include "constitutive/capillaryPressure/Layouts.hpp"
28 #include "constitutive/fluid/multifluid/MultiFluidBase.hpp"
29 #include "constitutive/solid/CoupledSolidBase.hpp"
30 #include "physicsSolvers/fluidFlow/kernels/compositional/AccumulationKernel.hpp"
33 
34 
35 namespace geos
36 {
37 
42 {
45 };
46 
51  "ComponentDensities",
52  "OverallComposition" );
53 
54 //START_SPHINX_INCLUDE_00
61 {
62 public:
63 
69  CompositionalMultiphaseBase( const string & name,
70  Group * const parent );
71 
74 
77 
80 
83 
86 
90  virtual ~CompositionalMultiphaseBase() override = default;
91 
92 //START_SPHINX_INCLUDE_01
93 
94  virtual void registerDataOnMesh( Group & meshBodies ) override;
95 
96  virtual void registerDataForCFL( Group & meshBodies ) { GEOS_UNUSED_VAR( meshBodies ); }
97 
105  virtual void
106  implicitStepSetup( real64 const & time_n,
107  real64 const & dt,
108  DomainPartition & domain ) override;
109 
110  virtual void
111  assembleSystem( real64 const time_n,
112  real64 const dt,
113  DomainPartition & domain,
114  DofManager const & dofManager,
115  CRSMatrixView< real64, globalIndex const > const & localMatrix,
116  arrayView1d< real64 > const & localRhs ) override;
117 
118  virtual void
120  real64 const dt,
121  DomainPartition & domain,
122  DofManager const & dofManager,
123  CRSMatrixView< real64, globalIndex const > const & localMatrix,
124  arrayView1d< real64 > const & localRhs ) override;
125 
135  template< typename SUBREGION_TYPE >
136  void accumulationAssemblyLaunch( DofManager const & dofManager,
137  SUBREGION_TYPE const & subRegion,
138  CRSMatrixView< real64, globalIndex const > const & localMatrix,
139  arrayView1d< real64 > const & localRhs ) const;
140 
141  virtual void
143 
144  virtual void
146  real64 const & dt,
147  DomainPartition & domain ) override;
148 
154 
160 
165  void updateFluidModel( ObjectManagerBase & dataGroup ) const;
166 
171  void updateRelPermModel( ObjectManagerBase & dataGroup ) const;
172 
177  void updateCapPressureModel( ObjectManagerBase & dataGroup ) const;
178 
183  void updateCompAmount( ElementSubRegionBase & subRegion ) const;
184 
189  void updateEnergy( ElementSubRegionBase & subRegion ) const;
190 
196 
201  virtual void updatePhaseMobility( ObjectManagerBase & dataGroup ) const = 0;
202 
203  real64 updateFluidState( ElementSubRegionBase & subRegion ) const;
204 
205  virtual void saveConvergedState( ElementSubRegionBase & subRegion ) const override final;
206 
207  virtual void saveSequentialIterationState( DomainPartition & domain ) override final;
208 
209  virtual void updateState( DomainPartition & domain ) override final;
210 
216 
221  integer numFluidPhases() const { return m_numPhases; }
222 
228 
232  virtual units::Unit getMassUnit() const override
233  { return m_useMass ? units::Unit::Mass : units::Unit::Mole; }
234 
245  DofManager const & dofManager,
246  CRSMatrixView< real64, globalIndex const > const & localMatrix,
247  arrayView1d< real64 > const & localRhs ) const;
248 
258  virtual void
260  DomainPartition const & domain,
261  DofManager const & dofManager,
262  CRSMatrixView< real64, globalIndex const > const & localMatrix,
263  arrayView1d< real64 > const & localRhs ) const = 0;
264 
273  virtual void
275  DomainPartition const & domain,
276  DofManager const & dofManager,
277  CRSMatrixView< real64, globalIndex const > const & localMatrix,
278  arrayView1d< real64 > const & localRhs ) const = 0;
279 
284  {
285  static constexpr char const * elemDofFieldString() { return "compositionalVariables"; }
286 
287  // inputs
288 
289  static constexpr char const * useMassFlagString() { return "useMass"; }
290  static constexpr char const * formulationTypeString() { return "formulationType"; }
291 
292  // time stepping controls
293 
294  static constexpr char const * solutionChangeScalingFactorString() { return "solutionChangeScalingFactor"; }
295  static constexpr char const * targetRelativePresChangeString() { return "targetRelativePressureChangeInTimeStep"; }
296  static constexpr char const * targetRelativeTempChangeString() { return "targetRelativeTemperatureChangeInTimeStep"; }
297  static constexpr char const * targetPhaseVolFracChangeString() { return "targetPhaseVolFractionChangeInTimeStep"; }
298  static constexpr char const * targetRelativeCompDensChangeString() { return "targetRelativeCompDensChangeInTimeStep"; }
299  static constexpr char const * targetCompFracChangeString() { return "targetCompFracChangeInTimeStep"; }
300  static constexpr char const * targetFlowCFLString() { return "targetFlowCFL"; }
301 
302 
303  // nonlinear solver parameters
304 
305  static constexpr char const * maxCompFracChangeString() { return "maxCompFractionChange"; }
306  static constexpr char const * maxRelativePresChangeString() { return "maxRelativePressureChange"; }
307  static constexpr char const * maxRelativeTempChangeString() { return "maxRelativeTemperatureChange"; }
308  static constexpr char const * maxRelativeCompDensChangeString() { return "maxRelativeCompDensChange"; }
309  static constexpr char const * allowLocalCompDensChoppingString() { return "allowLocalCompDensityChopping"; }
310  static constexpr char const * useTotalMassEquationString() { return "useTotalMassEquation"; }
311  static constexpr char const * useSimpleAccumulationString() { return "useSimpleAccumulation"; }
312  static constexpr char const * minCompDensString() { return "minCompDens"; }
313  static constexpr char const * minCompFracString() { return "minCompFrac"; }
314  static constexpr char const * maxSequentialCompDensChangeString() { return "maxSequentialCompDensChange"; }
315  static constexpr char const * minScalingFactorString() { return "minScalingFactor"; }
316 
317  static constexpr char const * relPermNamesString() { return "relPermNames"; }
318  static constexpr char const * capPressureNamesString() { return "capPressureNames"; }
319  static constexpr char const * diffusionNamesString() { return "diffusionNames"; }
320  static constexpr char const * dispersionNamesString() { return "dispersionNames"; }
321  };
322 
331  virtual void initializeFluidState( MeshLevel & mesh, string_array const & regionNames ) override;
332 
333  virtual void initializeThermalState( MeshLevel & mesh, string_array const & regionNames ) override;
334 
338  virtual void computeHydrostaticEquilibrium( DomainPartition & domain ) override;
339 
349  void applyDirichletBC( real64 const time,
350  real64 const dt,
351  DofManager const & dofManager,
352  DomainPartition & domain,
353  CRSMatrixView< real64, globalIndex const > const & localMatrix,
354  arrayView1d< real64 > const & localRhs ) const;
355 
365  void applySourceFluxBC( real64 const time,
366  real64 const dt,
367  DofManager const & dofManager,
368  DomainPartition & domain,
369  CRSMatrixView< real64, globalIndex const > const & localMatrix,
370  arrayView1d< real64 > const & localRhs ) const;
371 
381  virtual void applyAquiferBC( real64 const time,
382  real64 const dt,
383  DofManager const & dofManager,
384  DomainPartition & domain,
385  CRSMatrixView< real64, globalIndex const > const & localMatrix,
386  arrayView1d< real64 > const & localRhs ) const = 0;
387 
400  real64 const dt,
401  DofManager const & dofManager,
402  DomainPartition & domain,
403  CRSMatrixView< real64, globalIndex const > const & localMatrix,
404  arrayView1d< real64 > const & localRhs ) const;
405 
406 
412  void chopNegativeDensities( ElementSubRegionBase & subRegion );
413 
419 
420  virtual real64 setNextDtBasedOnStateChange( real64 const & currentDt,
421  DomainPartition & domain ) override;
422 
423 
424  virtual void computeCFLNumbers( DomainPartition & domain, real64 const & dt, real64 & maxPhaseCFL, real64 & maxCompCFL )
425  {
426  GEOS_UNUSED_VAR( domain, dt, maxPhaseCFL, maxCompCFL );
427  GEOS_ERROR( GEOS_FMT( "{}: computeCFLNumbers is not implemented for {}", getDataContext(), getCatalogName() ) );
428  }
429 
431 
432  integer useSimpleAccumulation() const { return m_useSimpleAccumulation; }
433 
434  integer useTotalMassEquation() const { return m_useTotalMassEquation; }
435 
436  virtual bool checkSequentialSolutionIncrements( DomainPartition & domain ) const override;
437 
438 protected:
439 
440  virtual void postInputInitialization() override;
441 
442  virtual void initializePreSubGroups() override;
443 
444 
451  void validateConstitutiveModels( DomainPartition const & domain ) const;
452 
457  void initializeAquiferBC( constitutive::ConstitutiveManager const & cm ) const;
458 
461 
464 
467 
470 
473 
476 
479 
482 
485 
488 
491 
494 
497 
500 
503 
506 
509 
512 
515 
518 
521 
524 
527 
530 
533  real64 m_maxSequentialCompDensChange;
534 
535 private:
536 
542  bool validateDirichletBC( DomainPartition & domain,
543  real64 const time ) const;
544 
545  virtual void setConstitutiveNames( ElementSubRegionBase & subRegion ) const override;
546 
547 };
548 
549 template< typename SUBREGION_TYPE >
551  SUBREGION_TYPE const & subRegion,
552  CRSMatrixView< real64, globalIndex const > const & localMatrix,
553  arrayView1d< real64 > const & localRhs ) const
554 {
555  constitutive::MultiFluidBase const & fluid =
556  getConstitutiveModel< constitutive::MultiFluidBase >( subRegion, subRegion.template getReference< string >( viewKeyStruct::fluidNamesString() ) );
557  constitutive::CoupledSolidBase const & solid =
558  getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, subRegion.template getReference< string >( viewKeyStruct::solidNamesString() ) );
559 
560  string const dofKey = dofManager.getKey( viewKeyStruct::elemDofFieldString() );
561 
562  using namespace isothermalCompositionalMultiphaseBaseKernels;
563 
564  BitFlags< KernelFlags > kernelFlags;
566  kernelFlags.set( KernelFlags::TotalMassEquation );
568  kernelFlags.set( KernelFlags::SimpleAccumulation );
569 
571  {
572  // isothermal for now
573  isothermalCompositionalMultiphaseBaseKernels::
574  AccumulationZFormulationKernelFactory::
575  createAndLaunch< parallelDevicePolicy<> >( m_numComponents,
576  m_numPhases,
577  dofManager.rankOffset(),
578  kernelFlags,
579  dofKey,
580  subRegion,
581  fluid,
582  solid,
583  localMatrix,
584  localRhs );
585  }
586  else
587  {
588  if( m_isThermal )
589  {
590  thermalCompositionalMultiphaseBaseKernels::
591  AccumulationKernelFactory::
592  createAndLaunch< parallelDevicePolicy<> >( m_numComponents,
593  m_numPhases,
594  dofManager.rankOffset(),
595  kernelFlags,
596  dofKey,
597  subRegion,
598  fluid,
599  solid,
600  localMatrix,
601  localRhs );
602  }
603  else
604  {
605  isothermalCompositionalMultiphaseBaseKernels::
606  AccumulationKernelFactory::
607  createAndLaunch< parallelDevicePolicy<> >( m_numComponents,
608  m_numPhases,
609  dofManager.rankOffset(),
610  kernelFlags,
611  dofKey,
612  subRegion,
613  fluid,
614  solid,
615  localMatrix,
616  localRhs );
617  }
618  }
619 }
620 
621 } // namespace geos
622 
623 #endif //GEOS_PHYSICSSOLVERS_FLUIDFLOW_COMPOSITIONALMULTIPHASEBASE_HPP_
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
#define GEOS_ERROR(...)
Raise a hard error and terminate the program.
Definition: Logger.hpp:226
Unit
Enumerator of available unit types for given physical scales. Units are in SI by default.
Definition: Units.hpp:59
void initializeAquiferBC(constitutive::ConstitutiveManager const &cm) const
Initialize the aquifer boundary condition (gravity vector, water phase index)
bool m_hasCapPressure
flag to determine whether or not to apply capillary pressure
real64 m_solutionChangeScalingFactor
damping factor for solution change targets
CompositionalMultiphaseBase(CompositionalMultiphaseBase &&)=default
default move constructor
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_targetCompFracChange
target (absolute) change in component fraction in a time step
virtual void initializeFluidState(MeshLevel &mesh, string_array const &regionNames) override
Initialize all variables from initial conditions.
real64 m_maxRelativePresChange
maximum (relative) change in pressure in a Newton iteration
real64 m_targetRelativePresChange
target (relative) change in pressure in a time step
real64 m_minCompFrac
minimum allowed global component fraction
void validateConstitutiveModels(DomainPartition const &domain) const
Utility function that checks the consistency of the constitutive models.
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.
bool m_hasDiffusion
flag to determine whether or not to apply diffusion
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
void chopNegativeCompFractions(DomainPartition &domain)
Sets all the negative component fractions (if any) to zero.
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
CompositionalMultiphaseFormulationType m_formulationType
formulation type
CompositionalMultiphaseBase(const string &name, Group *const parent)
main constructor for Group Objects
integer m_numPhases
the max number of fluid phases
bool m_hasDispersion
flag to determine whether or not to apply dispersion
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:45
globalIndex rankOffset(string const &fieldName) const
string const & getKey(string const &fieldName) const
Return the key used to record the field in the DofManager.
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
integer m_isThermal
flag to determine whether or not this is a thermal simulation
Class facilitating the representation of a multi-level discretization of a MeshBody.
Definition: MeshLevel.hpp:42
The ObjectManagerBase is the base object of all object managers in the mesh data hierachy.
virtual string getCatalogName() const =0
DataContext const & getDataContext() const
Definition: Group.hpp:1343
virtual void resetStateToBeginningOfStep(DomainPartition &domain) override
reset state of physics back to the beginning of the step.
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 assembleLocalTerms(DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const
assembles the accumulation other local terms for all cells
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.
void accumulationAssemblyLaunch(DofManager const &dofManager, SUBREGION_TYPE const &subRegion, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const
assembles the accumulation terms for all cells of a spcefici subRegion.
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:179
stdVector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:361
double real64
64-bit floating point type.
Definition: DataTypes.hpp:98
LvArray::CRSMatrixView< T, COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:309
CompositionalMultiphaseFormulationType
Options for flow formulation.
@ OverallComposition
use overall composition (z_c) as primary variables
@ ComponentDensities
use component densities as primary variables
int integer
Signed integer type.
Definition: DataTypes.hpp:81
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "richardson", "preconditioner")
Declare strings associated with enumeration values.