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 #include "constitutive/fluid/multifluid/MultiFluidBase.hpp"
26 #include "constitutive/solid/CoupledSolidBase.hpp"
27 #include "physicsSolvers/fluidFlow/kernels/compositional/AccumulationKernel.hpp"
29 
30 namespace geos
31 {
32 
37 {
40 };
41 
46  "ComponentDensities",
47  "OverallComposition" );
48 
49 //START_SPHINX_INCLUDE_00
56 {
57 public:
58 
64  CompositionalMultiphaseBase( const string & name,
65  Group * const parent );
66 
69 
72 
75 
78 
81 
85  virtual ~CompositionalMultiphaseBase() override = default;
86 
87 //START_SPHINX_INCLUDE_01
88 
89  virtual void registerDataOnMesh( Group & meshBodies ) override;
90 
91  virtual void registerDataForCFL( Group & meshBodies ) { GEOS_UNUSED_VAR( meshBodies ); }
92 
100  virtual void
101  implicitStepSetup( real64 const & time_n,
102  real64 const & dt,
103  DomainPartition & domain ) override;
104 
105  virtual void
106  assembleSystem( real64 const time_n,
107  real64 const dt,
108  DomainPartition & domain,
109  DofManager const & dofManager,
110  CRSMatrixView< real64, globalIndex const > const & localMatrix,
111  arrayView1d< real64 > const & localRhs ) override;
112 
113  virtual void
115  real64 const dt,
116  DomainPartition & domain,
117  DofManager const & dofManager,
118  CRSMatrixView< real64, globalIndex const > const & localMatrix,
119  arrayView1d< real64 > const & localRhs ) override;
120 
130  template< typename SUBREGION_TYPE >
131  void accumulationAssemblyLaunch( DofManager const & dofManager,
132  SUBREGION_TYPE const & subRegion,
133  CRSMatrixView< real64, globalIndex const > const & localMatrix,
134  arrayView1d< real64 > const & localRhs );
135 
136  virtual void
138 
139  virtual void
141  real64 const & dt,
142  DomainPartition & domain ) override;
143 
149 
155 
160  void updateFluidModel( ObjectManagerBase & dataGroup ) const;
161 
166  void updateRelPermModel( ObjectManagerBase & dataGroup ) const;
167 
172  void updateCapPressureModel( ObjectManagerBase & dataGroup ) const;
173 
178  void updateCompAmount( ElementSubRegionBase & subRegion ) const;
179 
184  void updateEnergy( ElementSubRegionBase & subRegion ) const;
185 
191 
196  virtual void updatePhaseMobility( ObjectManagerBase & dataGroup ) const = 0;
197 
198  real64 updateFluidState( ElementSubRegionBase & subRegion ) const;
199 
200  virtual void saveConvergedState( ElementSubRegionBase & subRegion ) const override final;
201 
202  virtual void saveSequentialIterationState( DomainPartition & domain ) override final;
203 
204  virtual void updateState( DomainPartition & domain ) override final;
205 
211 
216  integer numFluidPhases() const { return m_numPhases; }
217 
223 
227  virtual units::Unit getMassUnit() const override
228  { return m_useMass ? units::Unit::Mass : units::Unit::Mole; }
229 
240  DofManager const & dofManager,
241  CRSMatrixView< real64, globalIndex const > const & localMatrix,
242  arrayView1d< real64 > const & localRhs ) const;
243 
253  virtual void
255  DomainPartition const & domain,
256  DofManager const & dofManager,
257  CRSMatrixView< real64, globalIndex const > const & localMatrix,
258  arrayView1d< real64 > const & localRhs ) const = 0;
259 
268  virtual void
270  DomainPartition const & domain,
271  DofManager const & dofManager,
272  CRSMatrixView< real64, globalIndex const > const & localMatrix,
273  arrayView1d< real64 > const & localRhs ) const = 0;
274 
279  {
280  static constexpr char const * elemDofFieldString() { return "compositionalVariables"; }
281 
282  // inputs
283 
284  static constexpr char const * useMassFlagString() { return "useMass"; }
285  static constexpr char const * formulationTypeString() { return "formulationType"; }
286  static constexpr char const * relPermNamesString() { return "relPermNames"; }
287  static constexpr char const * capPressureNamesString() { return "capPressureNames"; }
288  static constexpr char const * diffusionNamesString() { return "diffusionNames"; }
289  static constexpr char const * dispersionNamesString() { return "dispersionNames"; }
290 
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  };
318 
327  virtual void initializeFluidState( MeshLevel & mesh, string_array const & regionNames ) override;
328 
329  virtual void initializeThermalState( MeshLevel & mesh, string_array const & regionNames ) override;
330 
334  virtual void computeHydrostaticEquilibrium( DomainPartition & domain ) override;
335 
345  void applyDirichletBC( real64 const time,
346  real64 const dt,
347  DofManager const & dofManager,
348  DomainPartition & domain,
349  CRSMatrixView< real64, globalIndex const > const & localMatrix,
350  arrayView1d< real64 > const & localRhs ) const;
351 
361  void applySourceFluxBC( real64 const time,
362  real64 const dt,
363  DofManager const & dofManager,
364  DomainPartition & domain,
365  CRSMatrixView< real64, globalIndex const > const & localMatrix,
366  arrayView1d< real64 > const & localRhs ) const;
367 
377  virtual void applyAquiferBC( real64 const time,
378  real64 const dt,
379  DofManager const & dofManager,
380  DomainPartition & domain,
381  CRSMatrixView< real64, globalIndex const > const & localMatrix,
382  arrayView1d< real64 > const & localRhs ) const = 0;
383 
396  real64 const dt,
397  DofManager const & dofManager,
398  DomainPartition & domain,
399  CRSMatrixView< real64, globalIndex const > const & localMatrix,
400  arrayView1d< real64 > const & localRhs ) const;
401 
402 
408  void chopNegativeDensities( ElementSubRegionBase & subRegion );
409 
415 
416  virtual real64 setNextDtBasedOnStateChange( real64 const & currentDt,
417  DomainPartition & domain ) override;
418 
419 
420  virtual void computeCFLNumbers( DomainPartition & domain, real64 const & dt, real64 & maxPhaseCFL, real64 & maxCompCFL )
421  {
422  GEOS_UNUSED_VAR( domain, dt, maxPhaseCFL, maxCompCFL );
423  GEOS_ERROR( GEOS_FMT( "{}: computeCFLNumbers is not implemented for {}", getDataContext(), getCatalogName() ) );
424  }
425 
427 
428  integer useSimpleAccumulation() const { return m_useSimpleAccumulation; }
429 
430  integer useTotalMassEquation() const { return m_useTotalMassEquation; }
431 
432  virtual bool checkSequentialSolutionIncrements( DomainPartition & domain ) const override;
433 
434 protected:
435 
436  virtual void postInputInitialization() override;
437 
438  virtual void initializePreSubGroups() override;
439 
440 
447  void validateConstitutiveModels( DomainPartition const & domain ) const;
448 
453  void initializeAquiferBC( constitutive::ConstitutiveManager const & cm ) const;
454 
464  template< typename OBJECT_TYPE >
465  void applyFieldValue( real64 const & time_n,
466  real64 const & dt,
467  MeshLevel & mesh,
468  char const logMessage[],
469  string const fieldKey,
470  string const boundaryFieldKey ) const;
471 
474 
477 
480 
483 
486 
489 
492 
495 
498 
501 
504 
507 
510 
513 
516 
519 
522 
525 
528 
531 
534 
537 
540 
543 
546  real64 m_maxSequentialCompDensChange;
547 
548 private:
549 
555  bool validateDirichletBC( DomainPartition & domain,
556  real64 const time ) const;
557 
558  virtual void setConstitutiveNames( ElementSubRegionBase & subRegion ) const override;
559 
560 };
561 
562 template< typename OBJECT_TYPE >
564  real64 const & dt,
565  MeshLevel & mesh,
566  char const logMessage[],
567  string const fieldKey,
568  string const boundaryFieldKey ) const
569 {
571 
572  fsManager.apply< OBJECT_TYPE >( time_n + dt,
573  mesh,
574  fieldKey,
575  [&]( FieldSpecificationBase const & fs,
576  string const & setName,
578  OBJECT_TYPE & targetGroup,
579  string const & )
580  {
581  if( fs.getLogLevel() >= 1 && m_nonlinearSolverParameters.m_numNewtonIterations == 0 )
582  {
583  globalIndex const numTargetElems = MpiWrapper::sum< globalIndex >( lset.size() );
584  GEOS_LOG_RANK_0( GEOS_FMT( logMessage,
585  getName(), time_n+dt, fs.getCatalogName(), fs.getName(),
586  setName, targetGroup.getName(), fs.getScale(), numTargetElems ) );
587  }
588 
589  // Specify the bc value of the field
590  fs.applyFieldValue< FieldSpecificationEqual,
591  parallelDevicePolicy<> >( lset,
592  time_n + dt,
593  targetGroup,
594  boundaryFieldKey );
595  } );
596 }
597 
598 template< typename SUBREGION_TYPE >
600  SUBREGION_TYPE const & subRegion,
601  CRSMatrixView< real64, globalIndex const > const & localMatrix,
602  arrayView1d< real64 > const & localRhs )
603 {
604  constitutive::MultiFluidBase const & fluid =
605  getConstitutiveModel< constitutive::MultiFluidBase >( subRegion, subRegion.template getReference< string >( viewKeyStruct::fluidNamesString() ) );
606  constitutive::CoupledSolidBase const & solid =
607  getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, subRegion.template getReference< string >( viewKeyStruct::solidNamesString() ) );
608 
609  string const dofKey = dofManager.getKey( viewKeyStruct::elemDofFieldString() );
610 
611  using namespace isothermalCompositionalMultiphaseBaseKernels;
612 
613  BitFlags< KernelFlags > kernelFlags;
615  kernelFlags.set( KernelFlags::TotalMassEquation );
617  kernelFlags.set( KernelFlags::SimpleAccumulation );
618 
619  if( m_isThermal )
620  {
621  thermalCompositionalMultiphaseBaseKernels::
622  AccumulationKernelFactory::
623  createAndLaunch< parallelDevicePolicy<> >( m_numComponents,
624  m_numPhases,
625  dofManager.rankOffset(),
626  kernelFlags,
627  dofKey,
628  subRegion,
629  fluid,
630  solid,
631  localMatrix,
632  localRhs );
633  }
634  else
635  {
636  isothermalCompositionalMultiphaseBaseKernels::
637  AccumulationKernelFactory::
638  createAndLaunch< parallelDevicePolicy<> >( m_numComponents,
639  m_numPhases,
640  dofManager.rankOffset(),
641  kernelFlags,
642  dofKey,
643  subRegion,
644  fluid,
645  solid,
646  localMatrix,
647  localRhs );
648  }
649 }
650 
651 } // namespace geos
652 
653 #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:56
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_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.
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
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
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
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...
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()
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
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:1345
string const & getName() const
Get group name.
Definition: Group.hpp:1331
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 accumulationAssemblyLaunch(DofManager const &dofManager, SUBREGION_TYPE const &subRegion, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)
assembles the accumulation terms for all cells of a spcefici subRegion.
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
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "preconditioner")
Declare strings associated with enumeration values.
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
CompositionalMultiphaseFormulationType
Options for flow formulation.
@ OverallComposition
use overall composition (z_c) as primary variables
@ ComponentDensities
use component densities as primary variables
std::vector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:393