pyscal.wateroil
Wateroil module
- class pyscal.wateroil.WaterOil(swirr=0.0, swl=0.0, swcr=0.0, sorw=0.0, socr=None, h=None, tag='', fast=False, _sgcr=None, _sgl=None)[source]
A representation of two-phase properties for oil-water.
Can hold relative permeability data, and capillary pressure.
- Parametrizations for relative permeability
Corey
LET
- For capillary pressure
Simplified J-function
For object initialization, only saturation endpoints must be inputted, and saturation resolution. An optional string can be added as a ‘tag’ that can be used when outputting.
Relative permeability and/or capillary pressure can be added through parametrizations, or from a dataframe (will incur interpolation).
Can be dumped as include files for Eclipse/OPM and Nexus simulators.
- Parameters:
swirr (
float
) – Absolute minimal water saturation at infinite capillary pressure.swl (
float
) – First water saturation point in generated table. Used for normalized saturations.swcr (
float
) – Critical water saturation. Water will not be mobile before the water saturation is above this value.sorw (
float
) – Residual oil saturation after water flooding. At this oil saturation, the oil has zero relative permeability.socr (
Optional
[float
]) – Critical oil saturation, oil will not be mobile before the oil saturation is above socr. This parameter will default to sorw and is to be used larger than sorw for oil paleo zone modelling cases.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
- SWFN(header=True, dataincommentrow=True, swcomment=None, crosspointcomment=None)[source]
Return a SWFN keyword with data to Eclipse
The columns sw, krw 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 SWFN string should be emitted. If you have multiple satnums, you should have True only for the first (or False for all, and emit the SWFN yourself). Defaults to True.dataincommentrow (
bool
) – boolean for wheter metadata should be printed, defaults to True.swcomment (
Optional
[str
]) – String to be used for swcomment, overrides what this object can provide. Used by GasWatercrosspointcomment (
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.
- Return type:
str
- SWOF(header=True, dataincommentrow=True)[source]
Produce SWOF input for Eclipse reservoir simulator.
The columns sw, krw, krow and pc are outputted and formatted accordingly.
Meta-information for the tabulated data are printed as Eclipse comments.
- Parameters:
header (
bool
) – Indicate whether the SWOF string should be emitted. If you have multiple SATNUMs, you should set this to True only for the first (or False for all, and emit the SWOF yourself). Default Truedataincommentrow (
bool
) – Wheter metadata should be printed. Defualt True
- Return type:
str
- WOTABLE(header=True, dataincommentrow=True)[source]
Return a string for a Nexus WOTABLE
- Return type:
str
- add_LET_oil(l=2.0, e=2.0, t=2.0, kroend=1)[source]
Add kro data through LET parametrization
- Parameters:
l (
float
) – LET parametere (
float
) – LET parametert (
float
) – LET parameterkroend (
float
) – value of kro at swl
- Return type:
None
- Returns:
None (modifies object)
- add_LET_pc_imb(Ls, Es, Ts, Lf, Ef, Tf, Pcmax, Pcmin, Pct)[source]
Add an imbition LET capillary pressure curve.
Docs: https://wiki.equinor.com/wiki/index.php/Res:The_LET_correlation_for_capillary_pressure
- Return type:
None
- add_LET_pc_pd(Lp, Ep, Tp, Lt, Et, Tt, Pcmax, Pct)[source]
Add a primary drainage LET capillary pressure curve.
Docs: https://wiki.equinor.com/wiki/index.php/Res:The_LET_correlation_for_capillary_pressure
Note that Pc where Sw > 1 - sorw will appear linear because there are no saturation points in that interval.
- Return type:
None
- add_LET_water(l=2.0, e=2.0, t=2.0, krwend=1.0, krwmax=None)[source]
Add krw data through LET parametrization
The LET model applies for sw < 1 - sgrw. For higher water saturations, krw is linear between krwend and krwmax.
krwmax will be ignored if sorw is close to zero.
- Parameters:
l (
float
) – LET parametere (
float
) – LET parametert (
float
) – LET parameterkrwend (
float
) – value of krw at 1 - sorwkrwmax (
Optional
[float
]) – maximal value at Sw=1. Default 1
- Return type:
None
- add_corey_oil(now=2.0, kroend=1.0)[source]
Add kro data through the Corey parametrization
Corey applies to the interval between swcr and 1 - sorw
Curve is linear between swl and swcr, zero above 1 - sorw.
- Parameters:
now (
float
) – Corey exponentkroend (
float
) – kro value at swcr
- Return type:
None
- Returns:
None (modifies object)
- add_corey_water(nw=2.0, krwend=1.0, krwmax=None)[source]
Add krw data through the Corey parametrization
A column named ‘krw’ will be added. If it exists, it will be replaced.
The Corey model applies for sw < 1 - sorw. For higher water saturations, krw is linear between krwend and krwmax.
krwmax will be ignored if sorw is close to zero
- Parameters:
nw (
float
) – Corey parameter for water.krwend (
float
) – value of krw at 1 - sorw.krwmax) – maximal value at Sw=1. Default 1
- Return type:
None
- add_fromtable(dframe, swcolname='SW', krwcolname='KRW', krowcolname='KROW', pccolname='PCOW', krwcomment='', krowcomment='', pccomment='', sorw=None, socr=None)[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 Sw-table
If you have krw and krow in different dataframes, call this function twice
Calling function is responsible for checking if any data was actually added to the table.
The relpermdata will be interpolated using a monotone cubic interpolator below 1-sorw, and linearly above 1-sorw. Capillary pressure data will be interpolated monotone cubicly over the entire saturation interval
The python package ecl2df has a tool for converting Eclipse input files to dataframes.
- Parameters:
dframe (
DataFrame
) – containing dataswcolname (
str
) – column name with the saturation data in the dataframe dfkrwcolname (
str
) – name of column in df with krwkrowcolname (
str
) – name of column in df with krowpccolname (
str
) – name of column in df with capillary pressure datakrwcomment (
str
) – Inserted into commentkrowcomment (
str
) – Inserted into commentpccomment (
str
) – Inserted into commentsorw (
Optional
[float
]) – Explicit sorw. If None, it will be estimated from the numbers in krwsocr (
Optional
[float
]) – Explicit socr. If None, it will be estimated from the numbers in krow
- Return type:
None
- add_normalized_J(a, b, poro, perm, sigma_costau)[source]
Add capillary pressure in bar through a normalized J-function.
\[p_c = \frac{\left(\frac{S_w}{a}\right)^{\frac{1}{b}} \sigma \cos \tau}{\sqrt{\frac{k}{\phi}}}\]The \(S_w\) saturation used in the formula is normalized with respect to the swirr parameter.
- Parameters:
a (
float
) – a parameterb (
float
) – b exponent (typically negative)poro (
float
) – Porosity value, fraction between 0 and 1perm (
float
) – Permeability value in mDsigma_costau (
float
) – Interfacial tension in mN/m (typical value 30 mN/m)
- Return type:
None
- Returns:
None. Modifies pc column in self.table, using bar as pressure unit.
- add_simple_J(a=5.0, b=-1.5, poro_ref=0.25, perm_ref=100, drho=300, g=9.81)[source]
Add capillary pressure function from a simplified J-function
This is the RMS version of the coefficients a and b, the formula used is
\[J = a S_w^b\]J is not dimensionless in this equation. The capillary pressure will be in bars.
This is identical to the also seen formula
\[J = 10^{b \log(S_w) + \log(a)}\]\(S_w\) in this formula is normalized with respect to the swirr variable of the WaterOil object.
- Parameters:
a (
float
) – a coefficientb (
float
) – b coefficientporo_ref (
float
) – Reference porosity for scaling to Pc, between 0 and 1perm_ref (
float
) – Reference permeability for scaling to Pc, in milliDarcydrho (
float
) – Density difference between water and oil, in SI units kg/m³. Default value is 300g (
float
) – Gravitational acceleration, in SI units m/s², default value is 9.81
- Return type:
None
- Returns:
None. Modifies pc column in self.table, using bar as pressure unit.
- add_simple_J_petro(a, b, poro_ref=0.25, perm_ref=100, drho=300, g=9.81)[source]
Add capillary pressure function from a simplified J-function
This is the petrophysical version of the coefficients a and b, the formula used is
\[J = \left(\frac{S_w}{a}\right)^{\frac{1}{b}}\]which is identical to
\[J = 10^\frac{\log(S_w) - \log(a)}{b}\]J is not dimensionless in this equation.
\(S_w\) in this formula is normalized with respect to the swirr variable of the WaterOil object.
- Parameters:
a (
float
) – a coefficient, petrophysical versionb (
float
) – b coefficient, petrophysical versionporo_ref (
float
) – Reference porosity for scaling to Pc, between 0 and 1perm_ref (
float
) – Reference permeability for scaling to Pc, in milliDarcydrho (
float
) – Density difference between water and oil, in SI units kg/m³. Default value is 300g (
float
) – Gravitational acceleration, in SI units m/s², default value is 9.81
- Return type:
None
- Returns:
None. Modifies pc column in self.table, using bar as pressure unit.
- add_skjaeveland_pc(cw, co, aw, ao, swr=None, sor=None)[source]
Add capillary pressure from the Skjæveland correlation,
Doc: https://wiki.equinor.com/wiki/index.php/Res:The_Skjaeveland_correlation_for_capillary_pressure
The implementation is unit independent, units are contained in the input constants.
If swr and sor are not provided, it will be taken from the swirr and sorw. Only use different values here if you know what you are doing.
Modifies or adds self.table[“PC”] if succesful. Returns false if error occured.
- crosspoint()[source]
Locate and return the saturation point where krw = krow
Accuracy of this crosspoint depends on the resolution chosen when initializing the saturation range
- Return type:
float
- Returns:
The water saturation where krw == krow, for relperm linearly interpolated in water saturation.
- estimate_sorw(curve='KRW')[source]
Estimate sorw of the current krw data.
This is mostly relevant when add_fromtable() has been used. sorw is estimated by searching for a linear part in krw downwards from sw=1. In practice it is impossible to infer sorw = 0, since we are limited by h, and the last segment from sw=1-h to sw=1 can always be assumed linear. Expect sorw = h if the real sorw = 0, but do not depend that it might not return zero in the future (one could argue that sorw = h should be specially treated to mean sorw = 0)
If the curve is linear everywhere, sorw will be returned as 1 - (swl + h)
krow is not used, and should probably not be, as it can be very close to zero before approaching sorw.
- Parameters:
curve (
str
) – Colum name of column to use, default is krw. If this is all linear, but krow is not, you might be better off with krow- Return type:
float
- Returns:
The estimated sorw.
- estimate_swcr(curve='KRW')[source]
Estimate swcr of the current krw data.
kwcr is estimated by searching for a linear part in krw upwards from sw=swl. In practice it is impossible to infer swcr = 0, since we are limited by h, and the first segment is assumed linear anyways.
If the curve is linear everywhere, swcr can end up at the right end of your saturation interval.
- Parameters:
curve (
str
) – Colum name of column to use, default is krow. If this is all linear, but krw is not, you might be better off with krw- Return type:
float
- Returns:
The estimated sgcr.
- plotkrwkrow(mpl_ax=None, color='blue', alpha=1, linewidth=1, linestyle='-', marker=None, label='', logyscale=False)[source]
Plot krw and krow
If the argument ‘mpl_ax’ is not supplied, a new plot window will be made. If supplied, it will draw on the specified axis.
- Return type:
None
- plotpc(mpl_ax=None, color='blue', alpha=1, linewidth=1, linestyle='-', label='', logyscale=False)[source]
Plot capillary pressure (pc)
If mpl_ax is supplied, the curve will be drawn on that, if not, a new axis (plot) will be made
- Return type:
None
- selfcheck(mode='SWOF')[source]
Check validities of the data in the table.
An unfinished object will return False.
If you call SWOF, this function must not return False
This function should not throw an exception, but capture the error and give an error.
- Parameters:
mode (
str
) – “SWOF” or “SWFN”. If SWFN, krow is not required.- Return type:
bool
- set_endpoints_linearpart_krow(kroend)[source]
Set linear parts of krow outside endpoints
Curve will be zero in [1 - socr, 1].
This function is used by add_corey_water(), and perhaps by other utility functions. It should not be necessary for end-users.
- Parameters:
kroend (
float
) – value of kro at swcr- Return type:
None
- set_endpoints_linearpart_krw(krwend, krwmax=None)[source]
Set linear parts of krw outside endpoints.
Curve will be linear from [1 - sorw, 1] (from krwmax to krwend) and zero in [swl, swcr].
This function is used by add_corey_water(), and perhaps by other utility functions. It should not be necessary for end-users.
- Parameters:
krwend (
float
) – krw at 1 - sorwrkrwmax (
Optional
[float
]) – krw at Sw=1. Default 1.
- Return type:
None