OpenCAEPoro  0.2.0 Sep/22/2022
A simulator for multicomponent porous media flow
Public Member Functions | Friends | List of all members
Reservoir Class Reference

#include <Reservoir.hpp>

Public Member Functions

void InputParam (ParamRead &param)
 
void Setup ()
 Setup static information for reservoir with input params.
 
void ApplyControl (const USI &i)
 Apply the control of ith critical time point.
 
void PrepareWell ()
 Calculate Well Properties at the beginning of each time step.
 
void CalWellFlux ()
 Calculate Flux between Bulk and Wells.
 
void CalWellTrans ()
 Calculate Trans of Wells.
 
void CalVpore ()
 Calculate pore of Bulks.
 
void CalKrPc ()
 Calculate Relative Permeability and Capillary for each Bulk.
 
void CalMaxChange ()
 Calculate Maximum Change of some reference variables for IMPEC.
 
void CalIPRT (const OCP_DBL &dt)
 Calculate num of Injection, Production.
 
OCP_INT CheckP (const bool &bulkCheck=true, const bool &wellCheck=true)
 Check if abnormal Pressure occurs.
 
bool CheckNi ()
 Check if abnormal Pressure occurs.
 
bool CheckVe (const OCP_DBL &Vlim) const
 Check error between Fluids and Pores.
 
OCP_USI GetBulkNum () const
 Return the num of Bulk.
 
OCP_USI GetMaxFIMBulk () const
 Return MaxNUMFIMBulk.
 
USI GetWellNum () const
 Return the num of Well.
 
USI GetComNum () const
 Return the num of Components.
 
void SetupWellBulk ()
 
void GetNTQT (const OCP_DBL &dt)
 
void AllocateAuxIMPEC ()
 Allocate memory for auxiliary variables used for IMPEC.
 
void InitIMPEC ()
 Initialize the properties of Reservoir for IMPEC.
 
OCP_DBL CalCFL (const OCP_DBL &dt)
 Calcluate the CFL number, including bulks and wells for IMPEC.
 
void CalFLuxIMPEC ()
 Calculate flux between bulks, bulks and wells.
 
void CalConnFluxIMPEC ()
 Calculate flux between bulks.
 
void MassConseveIMPEC (const OCP_DBL &dt)
 Calculate Ni according to Flux.
 
void CalFlashIMPEC ()
 Calculate Flash For IMPEC.
 
void UpdateLastStepIMPEC ()
 Update value of last step for IMPEC.
 
void AllocateMatIMPEC (LinearSystem &myLS) const
 Allocate Maxmimum memory for internal Matirx for IMPEC.
 
void AssembleMatIMPEC (LinearSystem &myLS, const OCP_DBL &dt) const
 Assemble Matrix for IMPEC.
 
void GetSolutionIMPEC (const vector< OCP_DBL > &u)
 Return the Solution to Reservoir Pressure for IMPEC.
 
void ResetWellIMPEC ()
 Reset Well for IMPEC.
 
void ResetVal01IMPEC ()
 Reset Capillary Pressure, Flux for IMPEC.
 
void ResetVal02IMPEC ()
 Reset Capillary Pressure, Moles of Componnets, Flux for IMPEC.
 
void ResetVal03IMPEC ()
 
void AllocateAuxFIM ()
 Allocate memory for auxiliary variables used for FIM.
 
void InitFIM ()
 Initialize the properties of Reservoir for FIM.
 
void InitFIM_n ()
 
void CalFlashDerivFIM ()
 Calculate Flash for FIM, some derivatives are needed.
 
void CalFlashDerivFIM_n ()
 
void CalKrPcDerivFIM ()
 Calculate Relative Permeability and Capillary and some derivatives for each Bulk.
 
void UpdateLastStepFIM ()
 Update value of last step for FIM.
 
void AllocateMatFIM (LinearSystem &myLS) const
 Allocate Maxmimum memory for internal Matirx for FIM.
 
void AssembleMatFIM (LinearSystem &myLS, const OCP_DBL &dt) const
 Assemble Matrix for FIM.
 
void AssembleMatFIM_n (LinearSystem &myLS, const OCP_DBL &dt) const
 
void GetSolutionFIM (const vector< OCP_DBL > &u, const OCP_DBL &dPmax, const OCP_DBL &dSmax)
 
void GetSolutionFIM_n (const vector< OCP_DBL > &u, const OCP_DBL &dPmax, const OCP_DBL &dSmax)
 
void GetSolution01FIM (const vector< OCP_DBL > &u)
 
void CalResFIM (ResFIM &resFIM, const OCP_DBL &dt)
 Calculate the Resiual for FIM, it's also RHS of Linear System.
 
void ResetFIM (const bool &flag)
 Reset FIM.
 
OCP_DBL GetNRdPmax ()
 Return NRdPmax.
 
OCP_DBL GetNRdSmax (OCP_USI &index)
 Return NRdSmax.
 
OCP_DBL GetNRdNmax ()
 Return NRdNmax.
 
OCP_DBL GetNRdSmaxP ()
 Return NRdSmaxP.
 
void PrintSolFIM (const string &outfile) const
 
void ShowRes (const vector< OCP_DBL > &res) const
 
void AllocateAuxAIMt ()
 Allocate memory for auxiliary variables used for AIMt.
 
void AllocateMatAIMt (LinearSystem &myLS) const
 Allocate Maxmimum memory for internal Matirx for local FIM.
 
void SetupFIMBulk (const bool &NRflag=false)
 Setup FIMBulk.
 
void AddFIMBulk ()
 
void SetupFIMBulkBoundAIMs ()
 
void CalFlashDerivAIM (const bool &IfAIMs)
 Calculate Flash for local FIM, some derivatives are needed.
 
void CalKrPcDerivAIM (const bool &IfAIMs)
 Calculate Relative Permeability and Capillary and some derivatives for each Bulk.
 
void CalResAIMt (ResFIM &resFIM, const OCP_DBL &dt)
 Calculate the Resiual for local FIM, it's also RHS of Linear System.
 
void AssembleMatAIMt (LinearSystem &myLS, const OCP_DBL &dt) const
 Assemble Matrix for AIMt -— local FIM here.
 
void GetSolutionAIMt (const vector< OCP_DBL > &u, const OCP_DBL &dPmax, const OCP_DBL &dSmax)
 
void AllocateAuxAIMs ()
 Allocate memory for auxiliary variables used for AIMs.
 
void CalResAIMs (ResFIM &resFIM, const OCP_DBL &dt)
 Calculate the Resiual for AIMs, it's also RHS of Linear System.
 
void AssembleMatAIMs (LinearSystem &myLS, vector< OCP_DBL > &res, const OCP_DBL &dt) const
 Assemble Matrix for AIMs.
 
void GetSolutionAIMs (const vector< OCP_DBL > &u, const OCP_DBL &dPmax, const OCP_DBL &dSmax)
 
void ResetValAIM ()
 
OCP_DBL CalCFLAIM (const OCP_DBL &dt)
 
void UpdateLastStepAIM ()
 Update value of last step for IMPEC.
 
void AllocateAuxAIMc ()
 Allocate memory for auxiliary variables used for FIM.
 
void AssembleMatAIMc (LinearSystem &myLS, const OCP_DBL &dt) const
 Assemble Matrix for AIMc.
 
void CalResAIMc (ResFIM &resFIM, const OCP_DBL &dt)
 Calculate the Resiual for FIM, it's also RHS of Linear System.
 
void CalFlashAIMc ()
 
void CalFlashAIMc01 ()
 
void CalKrPcAIMc ()
 
void CalFlashDerivAIMc ()
 Calculate Flash for local FIM, some derivatives are needed.
 
void CalKrPcDerivAIMc ()
 Calculate Relative Permeability and Capillary and some derivatives for each Bulk.
 
void GetSolutionAIMc (const vector< OCP_DBL > &u, const OCP_DBL &dPmax, const OCP_DBL &dSmax)
 
void InitAIMc ()
 
void UpdatePj ()
 

Friends

class OCPControl
 
class Summary
 
class CriticalInfo
 
class DetailInfo
 
class OCP_IMPEC
 
class OCP_FIM
 
class OCP_AIMs
 
class OCP_AIMt
 
class OCP_AIMc
 
class Solver
 

Detailed Description

Reservoir is the core component in our simulator, it contains the all reservoir information, and all operations on it.

Reservoir has four Core components. Grids contains the basic informations of all grids as a database of reservoir. Bulk only stores active grids, which defines the area used for calculation. AllWells contains the well information, it's used to manage operations related to wells. BulkConn contains connections between bulks(active grids).

Definition at line 30 of file Reservoir.hpp.

Member Function Documentation

◆ GetSolutionAIMt()

void Reservoir::GetSolutionAIMt ( const vector< OCP_DBL > &  u,
const OCP_DBL dPmax,
const OCP_DBL dSmax 
)

Return the Solution to Reservoir Pressure and moles of Components for FIM Exactly, it's a Newton step.

Definition at line 570 of file Reservoir.cpp.

572 {
573  bulk.GetSolAIMt(u, dPmax, dSmax);
574  allWells.GetSolAIMt(u, bulk.numFIMBulk, bulk.GetComNum() + 1);
575 }
void GetSolAIMt(const vector< OCP_DBL > &u, const OCP_USI &bId, const USI &len)
Get solution from solver class after linear system is solved for local FIM.
Definition: AllWells.cpp:692
void GetSolAIMt(const vector< OCP_DBL > &u, const OCP_DBL &dPmaxlim, const OCP_DBL &dSmaxlim)
Get the solution for local FIM after a Newton iteration.
Definition: Bulk.cpp:3360
USI GetComNum() const
Return the number of components.
Definition: Bulk.hpp:147

References Bulk::GetComNum(), Bulk::GetSolAIMt(), and AllWells::GetSolAIMt().

◆ GetSolutionFIM()

void Reservoir::GetSolutionFIM ( const vector< OCP_DBL > &  u,
const OCP_DBL dPmax,
const OCP_DBL dSmax 
)

Return the Solution to Reservoir Pressure and moles of Components for FIM Exactly, it's a Newton step.

Definition at line 409 of file Reservoir.cpp.

411 {
412  bulk.GetSolFIM(u, dPmax, dSmax);
413  allWells.GetSolFIM(u, bulk.GetBulkNum(), bulk.GetComNum() + 1);
414 }
void GetSolFIM(const vector< OCP_DBL > &u, const OCP_USI &bId, const USI &len)
Get solution from solver class after linear system is solved for FIM.
Definition: AllWells.cpp:545
OCP_USI GetBulkNum() const
Return the number of bulks.
Definition: Bulk.hpp:145
void GetSolFIM(const vector< OCP_DBL > &u, const OCP_DBL &dPmaxlim, const OCP_DBL &dSmaxlim)
Get the solution for FIM after a Newton iteration.
Definition: Bulk.cpp:2587

References Bulk::GetBulkNum(), Bulk::GetComNum(), Bulk::GetSolFIM(), and AllWells::GetSolFIM().

◆ InputParam()

void Reservoir::InputParam ( ParamRead param)

Input param from internal param data structure, which stores the params from input files.

Definition at line 18 of file Reservoir.cpp.

19 {
21 
22  grid.InputParam(param.paramRs);
23  bulk.InputParam(param.paramRs);
24  allWells.InputParam(param.paramWell);
25 
26 }
#define OCP_FUNCNAME
Print Function Name.
Definition: UtilError.hpp:73
void InputParam(const ParamWell &paramWell)
Input param from ParamWell.
Definition: AllWells.cpp:19
void InputParam(ParamReservoir &rs_param)
Input param from internal data structure ParamReservoir.
Definition: Bulk.cpp:24
void InputParam(const ParamReservoir &rs_param)
Input parameters from the internal param structure.
Definition: Grid.cpp:14
ParamWell paramWell
Read the well params.
Definition: ParamRead.hpp:43
ParamReservoir paramRs
Read the reservoir params.
Definition: ParamRead.hpp:42

References Grid::InputParam(), AllWells::InputParam(), Bulk::InputParam(), OCP_FUNCNAME, ParamRead::paramRs, and ParamRead::paramWell.

◆ ResetVal03IMPEC()

void Reservoir::ResetVal03IMPEC ( )

Reset Pressure, Capillary Pressure, Moles of components, Flux, Volume of Pores for IMPEC

Definition at line 275 of file Reservoir.cpp.

276 {
277  OCP_FUNCNAME;
278  bulk.ResetphaseNum();
279  bulk.ResetminEigenSkip();
280  bulk.ResetflagSkip();
281  bulk.ResetziSkip();
282  bulk.ResetPSkip();
283  bulk.ResetKs();
284 
285  bulk.ResetPj();
286  bulk.ResetNi();
287  bulk.ResetNt();
288  bulk.ResetFlash();
289  bulk.ResetVp();
290  conn.Reset();
291 
292  // Becareful! if recalculate the flash, result may be different because the initial
293  // flash was calculated by InitFlash not Flash.
294 }
void Reset()
Reset physcial values of the current step with the previous step.
Definition: BulkConn.cpp:137
void ResetPj()
Reset Pj to the ones of the last time step.
Definition: Bulk.hpp:199
void ResetNi()
Reset Ni to the ones of the last time step.
Definition: Bulk.hpp:201
void ResetVp()
Reset Vp to the ones of the last time step.
Definition: Bulk.hpp:205
void ResetNt()
Reset Nt to the ones of the last time step.
Definition: Bulk.hpp:194
void ResetFlash()
Reset variables in flash calculations.
Definition: Bulk.cpp:1880

References OCP_FUNCNAME, BulkConn::Reset(), Bulk::ResetFlash(), Bulk::ResetNi(), Bulk::ResetNt(), Bulk::ResetPj(), and Bulk::ResetVp().


The documentation for this class was generated from the following files: