GEOS
SinglePhaseBase.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_SINGLEPHASEBASE_HPP_
21 #define GEOS_PHYSICSSOLVERS_FLUIDFLOW_SINGLEPHASEBASE_HPP_
22 
26 #include "constitutive/fluid/singlefluid/SingleFluidBase.hpp"
27 #include "constitutive/solid/CoupledSolidBase.hpp"
28 
29 
30 namespace geos
31 {
32 
33 namespace constitutive
34 {
35 
36 class ConstitutiveBase;
37 
38 }
39 
46 {
47 public:
53  SinglePhaseBase( const string & name,
54  Group * const parent );
55 
56 
58  SinglePhaseBase() = delete;
59 
61  SinglePhaseBase( SinglePhaseBase const & ) = delete;
62 
65 
67  SinglePhaseBase & operator=( SinglePhaseBase const & ) = delete;
68 
71 
75  virtual ~SinglePhaseBase() override = default;
76 
77  virtual void registerDataOnMesh( Group & meshBodies ) override;
78 
85 
86  virtual void
87  implicitStepSetup( real64 const & time_n,
88  real64 const & dt,
89  DomainPartition & domain ) override;
90 
91  virtual void
92  assembleSystem( real64 const time_n,
93  real64 const dt,
94  DomainPartition & domain,
95  DofManager const & dofManager,
97  arrayView1d< real64 > const & localRhs ) override;
98 
99  virtual void
101  real64 const dt,
102  DomainPartition & domain,
103  DofManager const & dofManager,
104  CRSMatrixView< real64, globalIndex const > const & localMatrix,
105  arrayView1d< real64 > const & localRhs ) override;
106 
107  virtual real64
109  DofManager const & dofManager,
110  arrayView1d< real64 const > const & localSolution ) override;
111 
112  virtual bool
114  DofManager const & dofManager,
115  arrayView1d< real64 const > const & localSolution,
116  real64 const scalingFactor ) override;
117 
118  virtual void
120 
121  virtual void
123  real64 const & dt,
124  DomainPartition & domain ) override;
125 
127 
138  DofManager const & dofManager,
139  CRSMatrixView< real64, globalIndex const > const & localMatrix,
140  arrayView1d< real64 > const & localRhs );
141 
151  template< typename SUBREGION_TYPE >
152  void accumulationAssemblyLaunch( DofManager const & dofManager,
153  SUBREGION_TYPE const & subRegion,
154  CRSMatrixView< real64, globalIndex const > const & localMatrix,
155  arrayView1d< real64 > const & localRhs );
156 
166  virtual void
168  DomainPartition const & domain,
169  DofManager const & dofManager,
170  CRSMatrixView< real64, globalIndex const > const & localMatrix,
171  arrayView1d< real64 > const & localRhs ) = 0;
172 
182  virtual void
184  DomainPartition const & domain,
185  DofManager const & dofManager,
186  CRSMatrixView< real64, globalIndex const > const & localMatrix,
187  arrayView1d< real64 > const & localRhs ) = 0;
188 
199  virtual void
201  real64 const dt,
202  DomainPartition const & domain,
203  DofManager const & dofManager,
204  CRSMatrixView< real64, globalIndex const > const & localMatrix,
205  arrayView1d< real64 > const & localRhs,
206  string const & jumpDofKey ) = 0;
207 
218  virtual void
220  real64 const dt,
221  DomainPartition const & domain,
222  DofManager const & dofManager,
223  CRSMatrixView< real64, globalIndex const > const & localMatrix,
224  arrayView1d< real64 > const & localRhs,
225  CRSMatrixView< real64, localIndex const > const & dR_dAper ) = 0;
226 
228  {
229  static constexpr char const * elemDofFieldString() { return "singlePhaseVariables"; }
230  };
231 
241  void
242  applyDirichletBC( real64 const time_n,
243  real64 const dt,
244  DomainPartition & domain,
245  DofManager const & dofManager,
246  CRSMatrixView< real64, globalIndex const > const & localMatrix,
247  arrayView1d< real64 > const & localRhs ) const;
248 
258  void
259  applySourceFluxBC( real64 const time_n,
260  real64 const dt,
261  DomainPartition & domain,
262  DofManager const & dofManager,
263  CRSMatrixView< real64, globalIndex const > const & localMatrix,
264  arrayView1d< real64 > const & localRhs ) const;
265 
275  virtual void
276  applyAquiferBC( real64 const time,
277  real64 const dt,
278  DomainPartition & domain,
279  DofManager const & dofManager,
280  CRSMatrixView< real64, globalIndex const > const & localMatrix,
281  arrayView1d< real64 > const & localRhs ) const = 0;
282 
283  virtual void
284  updateState ( DomainPartition & domain ) override final;
285 
290  real64
292 
293 
298  virtual void
299  updateFluidModel( ObjectManagerBase & dataGroup ) const;
300 
305  void
306  updateMass( ElementSubRegionBase & subRegion ) const;
307 
312  void
313  updateEnergy( ElementSubRegionBase & subRegion ) const;
314 
320 
326 
331  void
332  updateMobility( ObjectManagerBase & dataGroup ) const;
333 
334  virtual void initializePreSubGroups() override;
335 
337 
338  virtual void initializeFluidState( MeshLevel & mesh, arrayView1d< string const > const & regionNames ) override;
339 
340  virtual void initializeThermalState( MeshLevel & mesh, arrayView1d< string const > const & regionNames ) override;
341 
345  virtual void computeHydrostaticEquilibrium( DomainPartition & domain ) override;
346 
350  virtual void updatePressureGradient( DomainPartition & domain )
351  {
352  GEOS_UNUSED_VAR( domain );
353  }
354 
367  real64 const dt,
368  DofManager const & dofManager,
369  DomainPartition & domain,
370  CRSMatrixView< real64, globalIndex const > const & localMatrix,
371  arrayView1d< real64 > const & localRhs ) const;
372 
373 protected:
374 
379  virtual void validateConstitutiveModels( DomainPartition & domain ) const;
380 
384  void initializeAquiferBC() const;
385 
386  virtual void setConstitutiveNamesCallSuper( ElementSubRegionBase & subRegion ) const override;
387 
392  virtual void saveConvergedState( ElementSubRegionBase & subRegion ) const override;
393 
398  {
405  };
406 
411  {
414  };
415 
426  virtual FluidPropViews getFluidProperties( constitutive::ConstitutiveBase const & fluid ) const;
427 
428  virtual ThermalFluidPropViews getThermalFluidProperties( constitutive::ConstitutiveBase const & fluid ) const;
429 
430 private:
431  virtual void setConstitutiveNames( ElementSubRegionBase & subRegion ) const override;
432 
433 };
434 
435 template< typename SUBREGION_TYPE >
437  SUBREGION_TYPE const & subRegion,
438  CRSMatrixView< real64, globalIndex const > const & localMatrix,
439  arrayView1d< real64 > const & localRhs )
440 {
441  geos::constitutive::SingleFluidBase const & fluid =
442  getConstitutiveModel< geos::constitutive::SingleFluidBase >( subRegion, subRegion.template getReference< string >( viewKeyStruct::fluidNamesString() ) );
443  //START_SPHINX_INCLUDE_COUPLEDSOLID
444  geos::constitutive::CoupledSolidBase const & solid =
445  getConstitutiveModel< geos::constitutive::CoupledSolidBase >( subRegion, subRegion.template getReference< string >( viewKeyStruct::solidNamesString() ) );
446  //END_SPHINX_INCLUDE_COUPLEDSOLID
447 
448  string const dofKey = dofManager.getKey( viewKeyStruct::elemDofFieldString() );
449 
450  if( m_isThermal )
451  {
452  thermalSinglePhaseBaseKernels::
453  AccumulationKernelFactory::
454  createAndLaunch< parallelDevicePolicy<> >( dofManager.rankOffset(),
455  dofKey,
456  subRegion,
457  fluid,
458  solid,
459  localMatrix,
460  localRhs );
461  }
462  else
463  {
464  singlePhaseBaseKernels::
465  AccumulationKernelFactory::
466  createAndLaunch< parallelDevicePolicy<> >( dofManager.rankOffset(),
467  dofKey,
468  subRegion,
469  fluid,
470  solid,
471  localMatrix,
472  localRhs );
473  }
474 }
475 
476 } /* namespace geos */
477 
478 #endif //GEOS_PHYSICSSOLVERS_FLUIDFLOW_SINGLEPHASEBASE_HPP_
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
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...
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 void computeHydrostaticEquilibrium(DomainPartition &domain) override
Compute the hydrostatic equilibrium using the compositions and temperature input tables.
virtual void updatePressureGradient(DomainPartition &domain)
Update the cell-wise pressure gradient.
virtual void setConstitutiveNamesCallSuper(ElementSubRegionBase &subRegion) const override
This function sets constitutive name fields on an ElementSubRegionBase, and calls the base function i...
SinglePhaseBase & operator=(SinglePhaseBase &&)=delete
deleted move operator
real64 updateFluidState(ElementSubRegionBase &subRegion) const
Function to update all constitutive state and dependent variables.
virtual void saveConvergedState(ElementSubRegionBase &subRegion) const override
Utility function to save the converged state.
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.
virtual ~SinglePhaseBase() override=default
default destructor
void updateEnergy(ElementSubRegionBase &subRegion) const
Function to update energy.
void initializeAquiferBC() const
Initialize the aquifer boundary condition (gravity vector, water phase index)
virtual void applyAquiferBC(real64 const time, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const =0
Apply aquifer boundary conditions to the system.
void updateMobility(ObjectManagerBase &dataGroup) const
Function to update fluid mobility.
virtual void assembleEDFMFluxTerms(real64 const time_n, real64 const dt, DomainPartition const &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs, string const &jumpDofKey)=0
assembles the flux terms for all cells for the poroelastic case
SinglePhaseBase & operator=(SinglePhaseBase const &)=delete
deleted assignment operator
void updateSolidInternalEnergyModel(ObjectManagerBase &dataGroup) const
Update all relevant solid internal energy models using current values of temperature.
SinglePhaseBase()=delete
deleted default constructor
SinglePhaseBase(SinglePhaseBase const &)=delete
deleted copy constructor
virtual void initializePostInitialConditionsPreSubGroups() override
Called by InitializePostInitialConditions() prior to initializing sub-Groups.
void updateMass(ElementSubRegionBase &subRegion) const
Function to update fluid mass.
void applySourceFluxBC(real64 const time_n, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const
Apply source flux boundary conditions to the system.
void applyDirichletBC(real64 const time_n, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) const
Function to perform the Application of Dirichlet type BC's.
virtual void validateConstitutiveModels(DomainPartition &domain) const
Checks constitutive models for consistency.
void assembleAccumulationTerms(DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)
assembles the accumulation terms for all cells
virtual void registerDataOnMesh(Group &meshBodies) override
Register wrappers that contain data on the mesh objects.
virtual void updateFluidModel(ObjectManagerBase &dataGroup) const
Function to update all constitutive models.
virtual void updateState(DomainPartition &domain) override final
Recompute all dependent quantities from primary variables (including constitutive models)
virtual FluidPropViews getFluidProperties(constitutive::ConstitutiveBase const &fluid) const
Extract properties from a fluid.
virtual void assembleStabilizedFluxTerms(real64 const dt, DomainPartition const &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)=0
assembles the flux terms for all cells including jump stabilization
virtual void assembleFluxTerms(real64 const dt, DomainPartition const &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)=0
assembles the flux terms for all cells
virtual void assembleHydrofracFluxTerms(real64 const time_n, real64 const dt, DomainPartition const &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs, CRSMatrixView< real64, localIndex const > const &dR_dAper)=0
assembles the flux terms for all cells for the hydrofracture case
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.
SinglePhaseBase(const string &name, Group *const parent)
main constructor for Group Objects
void updateThermalConductivity(ElementSubRegionBase &subRegion) const
Update thermal conductivity.
SinglePhaseBase(SinglePhaseBase &&)=default
default move constructor
virtual void initializePreSubGroups() override
Called by Initialize() prior to initializing sub-Groups.
virtual real64 scalingForSystemSolution(DomainPartition &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localSolution) override
Function to determine if the solution vector should be scaled back in order to maintain a known const...
virtual void implicitStepComplete(real64 const &time, real64 const &dt, DomainPartition &domain) override
perform cleanup for implicit timestep
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 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 bool checkSystemSolution(DomainPartition &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localSolution, real64 const scalingFactor) override
Function to check system solution for physical consistency and constraint violation.
virtual void implicitStepSetup(real64 const &time_n, real64 const &dt, DomainPartition &domain) override
function to perform setup for implicit timestep
virtual void resetStateToBeginningOfStep(DomainPartition &domain) override
reset state of physics back to the beginning of the step.
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
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:196
Structure holding views into fluid properties used by the base solver.
arrayView2d< real64 const > const dDens_dPres
derivative of density w.r.t. pressure
arrayView2d< real64 const > const visc
viscosity
real64 const defaultViscosity
default vi to use for new elements
real64 const defaultDensity
default density to use for new elements
arrayView2d< real64 const > const dVisc_dPres
derivative of viscosity w.r.t. pressure
arrayView2d< real64 const > const dens
density
Structure holding views into thermal fluid properties used by the base solver.
arrayView2d< real64 const > const dDens_dTemp
derivative of density w.r.t. temperature
arrayView2d< real64 const > const dVisc_dTemp
derivative of viscosity w.r.t. temperature