OpenCAEPoro
0.2.0 Sep/22/2022
A simulator for multicomponent porous media flow
|
#include <Reservoir.hpp>
Public Member Functions | |
void | InputParam (ParamRead ¶m) |
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 |
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.
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.
References Bulk::GetComNum(), Bulk::GetSolAIMt(), and AllWells::GetSolAIMt().
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.
References Bulk::GetBulkNum(), Bulk::GetComNum(), Bulk::GetSolFIM(), and AllWells::GetSolFIM().
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.
References Grid::InputParam(), AllWells::InputParam(), Bulk::InputParam(), OCP_FUNCNAME, ParamRead::paramRs, and ParamRead::paramWell.
void Reservoir::ResetVal03IMPEC | ( | ) |
Reset Pressure, Capillary Pressure, Moles of components, Flux, Volume of Pores for IMPEC
Definition at line 275 of file Reservoir.cpp.
References OCP_FUNCNAME, BulkConn::Reset(), Bulk::ResetFlash(), Bulk::ResetNi(), Bulk::ResetNt(), Bulk::ResetPj(), and Bulk::ResetVp().