pyscal.gasoil
Representing a GasOil object
- class pyscal.gasoil.GasOil(swirr=0.0, sgcr=0.0, h=None, swl=0.0, sorg=0.0, sgro=0.0, tag='', krgendanchor='sorg', fast=False, _sgl=None)[source]
Object to represent two-phase properties for gas and oil.
Parametrizations available for relative permeability:
Corey
LET
or data can alternatively be read in from tabulated data (as Pandas DataFrame).
No support (yet) to add capillary pressure.
krgend is by default anchored both to 1-swl-sorg, but can be set to anchor to 1-swl instead. If the krgendanchor argument is something else than the string sorg, it will be anchored to 1-swl.
- Parameters:
swirr (
float
) – Absolute minimal water saturation at infinite capillary pressure. Not in use currently, except for in informational headers and for consistency checks.swl (
float
) – First water saturation point in water tables. In GasOil, it is used to obtain the normalized oil and gas saturation.sgcr (
float
) – Critical gas saturation. Gas will not be mobile before the gas saturation is above this value.sorg (
float
) – Residual oil saturation after gas flooding. At this oil saturation, the oil has zero relative permeability.sgro (
float
) – Residual gas, for use in gas-condensate modelling. Used as an endpoint for parametrized oil curve. Must be zero or equal to sgcr for compatibility with Eclipse three-point scaling.krgendanchor (
str
) – Set to sorg (default) or something else, where to anchor krgend. If sorg, then the normalized gas saturation will be equal to 1 at 1 - swl - sorg, if not, it will be 1 at 1 - swl. If sorg is zero it does not matter. krgmax is only relevant when this anchor is sorg.h (
Optional
[float
]) – Saturation step-length in the outputted table.tag (
str
) – Optional string identifier, only used in comments.fast (
bool
) – Set to True if in order to skip some integrity checks and nice-to-have features. Not needed to set for normal pyscal runs, as speed is seldom crucial. Default False
- GOTABLE(header=True, dataincommentrow=True)[source]
Produce GOTABLE input for the Nexus reservoir simulator.
The columns sg, krg, krog and pc are outputted and formatted accordingly.
Meta-information for the tabulated data are printed as Eclipse comments.
- Parameters:
header (
bool
) – boolean for whether the SGOF string should be emitted. If you have multiple satnums, you should have True only for the first (or False for all, and emit the SGOF yourself). Defaults to True.dataincommentrow (
bool
) – boolean for wheter metadata should be printed, defaults to True.
- Return type:
str
- SGFN(header=True, dataincommentrow=True, sgcomment=None, crosspointcomment=None)[source]
Produce SGFN input for Eclipse reservoir simulator.
The columns sg, krg, and pc are outputted and formatted accordingly.
Meta-information for the tabulated data are printed as Eclipse comments.
- Parameters:
header (
bool
) – boolean for whether the SGFN string should be emitted. If you have multiple satnums, you should have True only for the first (or False for all, and emit the SGFN yourself). Defaults to True.dataincommentrow (
bool
) – boolean for wheter metadata should be printed, defaults to True.sgcomment (
Optional
[str
]) – Provide the string to include in the comment section for describing the saturation endpoints. Used by GasWater.crosspointcomment (
Optional
[str
]) – String to be used for crosspoint comment string, overrides what this object can provide. Used by GasWater. If None, it will be computed, use empty string to avoid.
- SGOF(header=True, dataincommentrow=True)[source]
Produce SGOF input for Eclipse reservoir simulator.
The columns sg, krg, krog and pc are outputted and formatted accordingly.
Meta-information for the tabulated data are printed as Eclipse comments.
- Parameters:
header (
bool
) – Whether the SGOF string should be emitted. If you have multiple satnums, you should have True only for the first (or False for all, and emit the SGOF yourself). Defaults to True.dataincommentrow (
bool
) – Whether metadata should be printed, defaults to True.
- Return type:
str
- SLGOF(header=True, dataincommentrow=True)[source]
Produce SLGOF input for Eclipse reservoir simulator.
The columns sl (liquid saturation), krg, krog and pc are outputted and formatted accordingly.
Meta-information for the tabulated data are printed as Eclipse comments.
- Parameters:
header (
bool
) – boolean for whether the SLGOF string should be emitted. If you have multiple satnums, you should have True only for the first (or False for all, and emit the SGOF yourself). Defaults to True.dataincommentrow (
bool
) – boolean for wheter metadata should be printed, defaults to True.
- Return type:
str
- add_LET_gas(l=2.0, e=2.0, t=2.0, krgend=1.0, krgmax=None)[source]
Add gas relative permability data through the LET parametrization
A column called ‘KRG’ will be added, replaced if it does not exist
If krgendanchor is sorg, the LET curve ends at krgend at sg = 1 - swl - sorg, and then linear up to krgmax at sg = 1 - swl. If not, it ends at krgend at sg = 1 - swl.
- Parameters:
l (
float
) – L parameter in LETe (
float
) – E parameter in LETt (
float
) – T parameter in LETkrgend (
float
) – Value of krg at normalized gas saturation 1krgmax (
Optional
[float
]) – Value of krg at gas saturation 1
- Returns:
None (modifies internal state)
- add_LET_oil(l=2.0, e=2.0, t=2.0, kroend=1.0, kromax=None)[source]
Add oil (vs gas) relative permeability data through the Corey parametrization.
A column named ‘krog’ will be added, replaced if it exists.
All values where sg > 1 - sorg - swl are set to zero.
- Parameters:
l (
float
) – L parametere (
float
) – E parametert (
float
) – T parameterkroend (
float
) – The value at gas saturation sgcrkromax (
Optional
[float
]) – Value at sg=0 for sgro > 0
- add_corey_gas(ng=2.0, krgend=1.0, krgmax=None)[source]
Add krg data through the Corey parametrization
A column called ‘krg’ will be added. If it exists, it will be replaced.
If krgendanchor is sorg, the Corey curve ends at krgend at sg = 1 - swl - sorg, and then linear up to krgmax at sg = 1 - swl. If not, it ends at krgend at sg = 1 - swl.
krgmax is only relevant if krgendanchor is ‘sorg’
- add_corey_oil(nog=2, kroend=1, kromax=None)[source]
Add kro data through the Corey parametrization
A column named ‘kro’ will be added to the internal DataFrame, replaced if it exists.
All values above 1 - sorg - swl are set to zero.
- Parameters:
nog (
float
) – Corey exponent for oilkroend (
float
) – Value for krog at normalized oil saturation 1kromax (
Optional
[float
]) – Value for sg=0 if sgro > 0.
- Returns:
None (modifies internal class state)
- add_fromtable(dframe, sgcolname='SG', krgcolname='KRG', krogcolname='KROG', pccolname='PCOG', krgcomment='', krogcomment='', pccomment='')[source]
Interpolate relpermdata from a dataframe.
The saturation range with endpoints must be set up beforehand, and must be compatible with the tabular input. The tabular input will be interpolated to the initialized Sg-table.
If you have krg and krog in different dataframes, call this function twice
Calling function is responsible for checking if any data was actually added to the table.
- crosspoint()[source]
Locate and return the saturation point where krg = krog
Accuracy of this crosspoint depends on the resolution chosen when initializing the saturation range (it uses linear interpolation to solve for the zero)
- Return type:
float
- Returns:
The gas saturation where krg == krog, for relperm linearly interpolated in gas saturation.
- estimate_sgcr()[source]
Estimate sgcr of the current krog data.
sgcr is the largest gas saturation for which the gas relative permeability is approximately zero, where approximately zero is roughly equivalent to how many digits are outputted by SGOF().
- Return type:
float
- Returns:
The estimated sgcr.
- estimate_sgro()[source]
Estimate sgro of the current krog data
sgro is estimated by searching for a linear part in kro from sg=0. In practice it is impossible to infer sgro = 0, since we are limited by h.
When initializing GasOil, sgro must be either zero or equal to sgcr. This function only estimates sgro, and will not guarantee that condition.
If the curve is linear everywhere, sgro will be returned as 1 - swl + h
- Returns:
The estimated sgro
- Return type:
float
- estimate_sorg()[source]
Estimate sorg of the current krg or krog data.
sorg is estimated by searching for a linear part in krg downwards from sg=1-swl. In practice it is impossible to infer sorg = 0, since we are limited by h, and the last segment from sg=1-swl-h to sg=1-swl must always be assumed linear.
If the curve is linear everywhere, sorg will be returned as sgcr + h
If krgend is anchored to sorg, krg data is used to infer sorg. If not, krg cannot be used for this, and krog is used. sorg might be overestimated when krog is used if it very close to zero before reaching sorg.
- Return type:
float
- Returns:
The estimated sorg.
- plotkrgkrog(mpl_ax=None, color='blue', alpha=1.0, linewidth=1, linestyle='-', marker=None, label=None, logyscale=False)[source]
Plot krg and krog
If mpl_ax is not None, it will be used as a matplotlib axis to plot on, if None, a fresh plot will be made.
- selfcheck(mode='SGOF')[source]
Check validities of the data in the table.
This is to catch errors that are either physically wrong or at least causes Eclipse 100 to stop.
Returns True if no errors are found, False if at least one is found.
If you call SGOF/SLGOF, this function must not return False.
- Parameters:
mode (
str
) – If mode is “SGFN”, krog is not required.- Return type:
bool
- set_endpoints_linearpart_krg(krgend, krgmax=None)[source]
Set linear parts of krg outside endpoints.
Curve is set to zero in [0, sgcr].
Given the default krgendanchor==sorg, the curve will be linear in [1 - swl - sorg, 1 - swl] (from krgend to krgmax). If not anchored to sorg, there is no linear part near sg=1-swl.
If krgendanchor is set to sorg (default), then the normalized gas saturation sgn (which is what is raised to the power of ng) is 1 at sg = 1 - swl - sorg. If not, it is 1 at sg = 1 - swl.
krgmax is only relevant if krgendanchor is ‘sorg’
This function is used by add_corey/LET_gas(), and perhaps by other utility functions. It should not be necessary for end-users.
- Parameters:
krgend (
float
) – krg at sg = 1 - swl - sorg.krgmax (
Optional
[float
]) – krg at Sg = 1 - swl. Default 1.
- set_endpoints_linearpart_krog(kroend, kromax=None)[source]
Set linear parts of krog outside endpoints.
Linear for sg in [0, sgro], from kromax to kroend, but nonzero sgro should only be used in gas-condensate modelling. When sgro is zero, kromax will be ignored.
Zero for sg above 1 - sorg - swl.
This function is used by add_corey/LET_oil(), and perhaps by other utility functions. It should not be necessary for end-users.
- Parameters:
kroend (
float
) – krog at sg=sgro, also if sgro=0kromax (
Optional
[float
]) – krog at sg=0 for sgro > 0