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/ConstitutiveBase.hpp"
27 #include "constitutive/fluid/singlefluid/SingleFluidLayouts.hpp"
28 #include "constitutive/fluid/singlefluid/SingleFluidUtils.hpp"
29 
30 
31 namespace geos
32 {
33 
34 namespace constitutive
35 {
36 
37 class ConstitutiveBase;
38 
39 }
40 
47 {
48 public:
49  using SingleFluidProp = constitutive::SingleFluidVar< real64, 2, constitutive::singlefluid::LAYOUT_FLUID, constitutive::singlefluid::LAYOUT_FLUID_DER >;
50 
56  SinglePhaseBase( const string & name,
57  Group * const parent );
58 
59 
61  SinglePhaseBase() = delete;
62 
64  SinglePhaseBase( SinglePhaseBase const & ) = delete;
65 
68 
70  SinglePhaseBase & operator=( SinglePhaseBase const & ) = delete;
71 
74 
78  virtual ~SinglePhaseBase() override
79  {}
80 
81  virtual void registerDataOnMesh( Group & meshBodies ) override;
82 
89 
90  virtual void
91  implicitStepSetup( real64 const & time_n,
92  real64 const & dt,
93  DomainPartition & domain ) override;
94 
95  virtual void
96  assembleSystem( real64 const time_n,
97  real64 const dt,
98  DomainPartition & domain,
99  DofManager const & dofManager,
100  CRSMatrixView< real64, globalIndex const > const & localMatrix,
101  arrayView1d< real64 > const & localRhs ) override;
102 
103  virtual void
105  real64 const dt,
106  DomainPartition & domain,
107  DofManager const & dofManager,
108  CRSMatrixView< real64, globalIndex const > const & localMatrix,
109  arrayView1d< real64 > const & localRhs ) override;
110 
111  virtual real64
113  DofManager const & dofManager,
114  arrayView1d< real64 const > const & localSolution ) override;
115 
116  virtual bool
118  DofManager const & dofManager,
119  arrayView1d< real64 const > const & localSolution,
120  real64 const scalingFactor ) override;
121 
122  virtual void
124 
125  virtual void
127  real64 const & dt,
128  DomainPartition & domain ) override;
129 
131 
142  DofManager const & dofManager,
143  CRSMatrixView< real64, globalIndex const > const & localMatrix,
144  arrayView1d< real64 > const & localRhs );
145 
155  template< typename SUBREGION_TYPE >
156  void accumulationAssemblyLaunch( DofManager const & dofManager,
157  SUBREGION_TYPE const & subRegion,
158  CRSMatrixView< real64, globalIndex const > const & localMatrix,
159  arrayView1d< real64 > const & localRhs );
160 
170  virtual void
172  DomainPartition const & domain,
173  DofManager const & dofManager,
174  CRSMatrixView< real64, globalIndex const > const & localMatrix,
175  arrayView1d< real64 > const & localRhs ) = 0;
176 
186  virtual void
188  DomainPartition const & domain,
189  DofManager const & dofManager,
190  CRSMatrixView< real64, globalIndex const > const & localMatrix,
191  arrayView1d< real64 > const & localRhs ) = 0;
192 
203  virtual void
205  real64 const dt,
206  DomainPartition const & domain,
207  DofManager const & dofManager,
208  CRSMatrixView< real64, globalIndex const > const & localMatrix,
209  arrayView1d< real64 > const & localRhs,
210  string const & jumpDofKey ) = 0;
211 
213  {
214  static constexpr char const * elemDofFieldString() { return "singlePhaseVariables"; }
215  };
216 
226  void
227  applyDirichletBC( real64 const time_n,
228  real64 const dt,
229  DomainPartition & domain,
230  DofManager const & dofManager,
231  CRSMatrixView< real64, globalIndex const > const & localMatrix,
232  arrayView1d< real64 > const & localRhs ) const;
233 
243  void
244  applySourceFluxBC( real64 const time_n,
245  real64 const dt,
246  DomainPartition & domain,
247  DofManager const & dofManager,
248  CRSMatrixView< real64, globalIndex const > const & localMatrix,
249  arrayView1d< real64 > const & localRhs ) const;
250 
260  virtual void
261  applyAquiferBC( real64 const time,
262  real64 const dt,
263  DomainPartition & domain,
264  DofManager const & dofManager,
265  CRSMatrixView< real64, globalIndex const > const & localMatrix,
266  arrayView1d< real64 > const & localRhs ) const = 0;
267 
268  virtual void
269  updateState ( DomainPartition & domain ) override final;
270 
275  real64
277 
278 
283  virtual void
284  updateFluidModel( ObjectManagerBase & dataGroup ) const;
285 
290  void
291  updateMass( ElementSubRegionBase & subRegion ) const;
292 
297  template< integer IS_THERMAL >
298  void updateMass( ElementSubRegionBase & subRegion ) const;
299 
304  void
305  updateEnergy( ElementSubRegionBase & subRegion ) const;
306 
312 
318 
323  void
324  updateMobility( ObjectManagerBase & dataGroup ) const;
325 
326  virtual void initializePreSubGroups() override;
327 
329 
330  virtual void initializeFluidState( MeshLevel & mesh, string_array const & regionNames ) override;
331 
332  virtual void initializeThermalState( MeshLevel & mesh, string_array const & regionNames ) override;
333 
337  virtual void computeHydrostaticEquilibrium( DomainPartition & domain ) override;
338 
342  virtual void updatePressureGradient( DomainPartition & domain )
343  {
344  GEOS_UNUSED_VAR( domain );
345  }
346 
359  real64 const dt,
360  DofManager const & dofManager,
361  DomainPartition & domain,
362  CRSMatrixView< real64, globalIndex const > const & localMatrix,
363  arrayView1d< real64 > const & localRhs ) const;
364 
365  void applyDeltaVolume( ElementSubRegionBase & subRegion ) const;
366 
367 protected:
368 
373  virtual void validateConstitutiveModels( DomainPartition & domain ) const;
374 
378  void initializeAquiferBC() const;
379 
380  virtual void setConstitutiveNamesCallSuper( ElementSubRegionBase & subRegion ) const override;
381 
386  virtual void saveConvergedState( ElementSubRegionBase & subRegion ) const override;
387 
392  {
399  };
400 
401 
412  virtual FluidPropViews getFluidProperties( constitutive::ConstitutiveBase const & fluid ) const;
413 
414 
415 private:
416  virtual void setConstitutiveNames( ElementSubRegionBase & subRegion ) const override;
417 
418 };
419 
420 template< typename SUBREGION_TYPE >
422  SUBREGION_TYPE const & subRegion,
423  CRSMatrixView< real64, globalIndex const > const & localMatrix,
424  arrayView1d< real64 > const & localRhs )
425 {
426  string const dofKey = dofManager.getKey( viewKeyStruct::elemDofFieldString() );
427 
428  if( m_isThermal )
429  {
430  thermalSinglePhaseBaseKernels::
431  AccumulationKernelFactory::
432  createAndLaunch< parallelDevicePolicy<> >( dofManager.rankOffset(),
433  dofKey,
434  subRegion,
435  localMatrix,
436  localRhs );
437  }
438  else
439  {
440  singlePhaseBaseKernels::
441  AccumulationKernelFactory::
442  createAndLaunch< parallelDevicePolicy<> >( dofManager.rankOffset(),
443  dofKey,
444  subRegion,
445  localMatrix,
446  localRhs );
447  }
448 }
449 
450 } /* namespace geos */
451 
452 #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.
void updateMass(ElementSubRegionBase &subRegion) const
Template function to update fluid mass.
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.
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 ~SinglePhaseBase() override
default destructor
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
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:188
stdVector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:401
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:318
double real64
64-bit floating point type.
Definition: DataTypes.hpp:98
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:204
ArrayView< T, 3, USD > arrayView3d
Alias for 3D array view.
Definition: DataTypes.hpp:220
Structure holding views into fluid properties used by the base solver.
arrayView3d< real64 const, constitutive::singlefluid::USD_FLUID_DER > const dDens
density derivatives
real64 const defaultViscosity
default vi to use for new elements
real64 const defaultDensity
default density to use for new elements
arrayView2d< real64 const, constitutive::singlefluid::USD_FLUID > const visc
viscosity
arrayView2d< real64 const, constitutive::singlefluid::USD_FLUID > const dens
density
arrayView3d< real64 const, constitutive::singlefluid::USD_FLUID_DER > const dVisc
viscosity derivatives