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 = default;
79 
80  virtual void registerDataOnMesh( Group & meshBodies ) override;
81 
88 
89  virtual void
90  implicitStepSetup( real64 const & time_n,
91  real64 const & dt,
92  DomainPartition & domain ) override;
93 
94  virtual void
95  assembleSystem( real64 const time_n,
96  real64 const dt,
97  DomainPartition & domain,
98  DofManager const & dofManager,
100  arrayView1d< real64 > const & localRhs ) override;
101 
102  virtual void
104  real64 const dt,
105  DomainPartition & domain,
106  DofManager const & dofManager,
107  CRSMatrixView< real64, globalIndex const > const & localMatrix,
108  arrayView1d< real64 > const & localRhs ) override;
109 
110  virtual real64
112  DofManager const & dofManager,
113  arrayView1d< real64 const > const & localSolution ) override;
114 
115  virtual bool
117  DofManager const & dofManager,
118  arrayView1d< real64 const > const & localSolution,
119  real64 const scalingFactor ) override;
120 
121  virtual void
123 
124  virtual void
126  real64 const & dt,
127  DomainPartition & domain ) override;
128 
130 
141  DofManager const & dofManager,
142  CRSMatrixView< real64, globalIndex const > const & localMatrix,
143  arrayView1d< real64 > const & localRhs );
144 
154  template< typename SUBREGION_TYPE >
155  void accumulationAssemblyLaunch( DofManager const & dofManager,
156  SUBREGION_TYPE const & subRegion,
157  CRSMatrixView< real64, globalIndex const > const & localMatrix,
158  arrayView1d< real64 > const & localRhs );
159 
169  virtual void
171  DomainPartition const & domain,
172  DofManager const & dofManager,
173  CRSMatrixView< real64, globalIndex const > const & localMatrix,
174  arrayView1d< real64 > const & localRhs ) = 0;
175 
185  virtual void
187  DomainPartition const & domain,
188  DofManager const & dofManager,
189  CRSMatrixView< real64, globalIndex const > const & localMatrix,
190  arrayView1d< real64 > const & localRhs ) = 0;
191 
202  virtual void
204  real64 const dt,
205  DomainPartition const & domain,
206  DofManager const & dofManager,
207  CRSMatrixView< real64, globalIndex const > const & localMatrix,
208  arrayView1d< real64 > const & localRhs,
209  string const & jumpDofKey ) = 0;
210 
212  {
213  static constexpr char const * elemDofFieldString() { return "singlePhaseVariables"; }
214  };
215 
225  void
226  applyDirichletBC( real64 const time_n,
227  real64 const dt,
228  DomainPartition & domain,
229  DofManager const & dofManager,
230  CRSMatrixView< real64, globalIndex const > const & localMatrix,
231  arrayView1d< real64 > const & localRhs ) const;
232 
242  void
243  applySourceFluxBC( real64 const time_n,
244  real64 const dt,
245  DomainPartition & domain,
246  DofManager const & dofManager,
247  CRSMatrixView< real64, globalIndex const > const & localMatrix,
248  arrayView1d< real64 > const & localRhs ) const;
249 
259  virtual void
260  applyAquiferBC( real64 const time,
261  real64 const dt,
262  DomainPartition & domain,
263  DofManager const & dofManager,
264  CRSMatrixView< real64, globalIndex const > const & localMatrix,
265  arrayView1d< real64 > const & localRhs ) const = 0;
266 
267  virtual void
268  updateState ( DomainPartition & domain ) override final;
269 
274  real64
276 
277 
282  virtual void
283  updateFluidModel( ObjectManagerBase & dataGroup ) const;
284 
289  void
290  updateMass( ElementSubRegionBase & subRegion ) const;
291 
296  template< integer IS_THERMAL >
297  void updateMass( ElementSubRegionBase & subRegion ) const;
298 
303  void
304  updateEnergy( ElementSubRegionBase & subRegion ) const;
305 
311 
317 
322  void
323  updateMobility( ObjectManagerBase & dataGroup ) const;
324 
325  virtual void initializePreSubGroups() override;
326 
328 
329  virtual void initializeFluidState( MeshLevel & mesh, string_array const & regionNames ) override;
330 
331  virtual void initializeThermalState( MeshLevel & mesh, string_array const & regionNames ) override;
332 
336  virtual void computeHydrostaticEquilibrium( DomainPartition & domain ) override;
337 
341  virtual void updatePressureGradient( DomainPartition & domain )
342  {
343  GEOS_UNUSED_VAR( domain );
344  }
345 
358  real64 const dt,
359  DofManager const & dofManager,
360  DomainPartition & domain,
361  CRSMatrixView< real64, globalIndex const > const & localMatrix,
362  arrayView1d< real64 > const & localRhs ) const;
363 
364  void applyDeltaVolume( ElementSubRegionBase & subRegion ) const;
365 
366 protected:
367 
372  virtual void validateConstitutiveModels( DomainPartition & domain ) const;
373 
377  void initializeAquiferBC() const;
378 
379  virtual void setConstitutiveNamesCallSuper( ElementSubRegionBase & subRegion ) const override;
380 
385  virtual void saveConvergedState( ElementSubRegionBase & subRegion ) const override;
386 
391  {
398  };
399 
400 
411  virtual FluidPropViews getFluidProperties( constitutive::ConstitutiveBase const & fluid ) const;
412 
413 
414 private:
415  virtual void setConstitutiveNames( ElementSubRegionBase & subRegion ) const override;
416 
417 };
418 
419 template< typename SUBREGION_TYPE >
421  SUBREGION_TYPE const & subRegion,
422  CRSMatrixView< real64, globalIndex const > const & localMatrix,
423  arrayView1d< real64 > const & localRhs )
424 {
425  string const dofKey = dofManager.getKey( viewKeyStruct::elemDofFieldString() );
426 
427  if( m_isThermal )
428  {
429  thermalSinglePhaseBaseKernels::
430  AccumulationKernelFactory::
431  createAndLaunch< parallelDevicePolicy<> >( dofManager.rankOffset(),
432  dofKey,
433  subRegion,
434  localMatrix,
435  localRhs );
436  }
437  else
438  {
439  singlePhaseBaseKernels::
440  AccumulationKernelFactory::
441  createAndLaunch< parallelDevicePolicy<> >( dofManager.rankOffset(),
442  dofKey,
443  subRegion,
444  localMatrix,
445  localRhs );
446  }
447 }
448 
449 } /* namespace geos */
450 
451 #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.
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
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
std::vector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:393
ArrayView< T, 3, USD > arrayView3d
Alias for 3D array view.
Definition: DataTypes.hpp:212
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