gstools.covmodel.base¶
GStools subpackage providing the base class for covariance models.
The following classes are provided
CovModel ([dim, var, len_scale, nugget, …]) |
Base class for the GSTools covariance models. |
-
class
gstools.covmodel.base.
CovModel
(dim=3, var=1.0, len_scale=1.0, nugget=0.0, anis=1.0, angles=0.0, integral_scale=None, var_raw=None, hankel_kw=None, **opt_arg)[source]¶ Bases:
object
Base class for the GSTools covariance models.
Parameters: - dim (
int
, optional) – dimension of the model. Default:3
- var (
float
, optional) – variance of the model (the nugget is not included in “this” variance) Default:1.0
- len_scale (
float
orlist
, optional) – length scale of the model. If a single value is given, the same length-scale will be used for every direction. If multiple values (for main and transversal directions) are given, anis will be recalculated accordingly. Default:1.0
- nugget (
float
, optional) – nugget of the model. Default:0.0
- anis (
float
orlist
, optional) – anisotropy ratios in the transversal directions [y, z]. Default:1.0
- angles (
float
orlist
, optional) –angles of rotation:
- in 2D: given as rotation around z-axis
- in 3D: given by yaw, pitch, and roll (known as Tait–Bryan angles)
Default:
0.0
- integral_scale (
float
orlist
orNone
, optional) – If given,len_scale
will be ignored and recalculated, so that the integral scale of the model matches the given one. Default:None
- var_raw (
float
orNone
, optional) – raw variance of the model which will be multiplied withCovModel.var_factor
to result in the actual variance. If given,var
will be ignored. (This is just for models that overrideCovModel.var_factor
) Default:None
- hankel_kw (
dict
orNone
, optional) – Modify the init-arguments ofhankel.SymmetricFourierTransform
used for the spectrum calculation. Use with caution (Better: Don’t!).None
is equivalent to{"a": -1, "b": 1, "N": 1000, "h": 0.001}
. Default:None
Examples
>>> from gstools import CovModel >>> import numpy as np >>> class Gau(CovModel): ... def cor(self, h): ... return np.exp(-h**2) ... >>> model = Gau() >>> model.spectrum(2) 0.00825830126008459
Attributes: angles
numpy.ndarray
: Rotation angles (in rad) of the model.anis
numpy.ndarray
: The anisotropy factors of the model.arg
arg_bounds
dict
: Bounds for all parameters.dim
int
: The dimension of the model.dist_func
do_rotation
bool
: State if a rotation is performed.hankel_kw
dict
:hankel.SymmetricFourierTransform
kwargs.has_cdf
bool
: State if a cdf is defined by the user.has_ppf
bool
: State if a ppf is defined by the user.integral_scale
float
: The main integral scale of the model.integral_scale_vec
numpy.ndarray
: The integral scales in each direction.len_scale
float
: The main length scale of the model.len_scale_bounds
list
: Bounds for the lenght scale.len_scale_vec
numpy.ndarray
: The length scales in each direction.name
str
: The name of the CovModel class.nugget
float
: The nugget of the model.nugget_bounds
list
: Bounds for the nugget.opt_arg
opt_arg_bounds
dict
: Bounds for the optional arguments.pykrige_angle
2D rotation angle for pykrige.
pykrige_angle_x
3D rotation angle around x for pykrige.
pykrige_angle_y
3D rotation angle around y for pykrige.
pykrige_angle_z
3D rotation angle around z for pykrige.
pykrige_anis
2D anisotropy ratio for pykrige.
pykrige_anis_y
3D anisotropy ratio in y direction for pykrige.
pykrige_anis_z
3D anisotropy ratio in z direction for pykrige.
pykrige_kwargs
Keyword arguments for pykrige routines.
sill
float
: The sill of the variogram.var
float
: The variance of the model.var_bounds
list
: Bounds for the variance.var_raw
float
: The raw variance of the model without factor.
Methods
calc_integral_scale
()Calculate the integral scale of the isotrope model. check_arg_bounds
()Check arguments to be within the given bounds. check_opt_arg
()Run checks for the optional arguments. cor_spatial
(pos)Spatial correlation respecting anisotropy and rotation. cov_nugget
(r)Covariance of the model respecting the nugget at r=0. cov_spatial
(pos)Spatial covariance respecting anisotropy and rotation. default_arg_bounds
()Provide default boundaries for arguments. default_opt_arg
()Provide default optional arguments by the user. default_opt_arg_bounds
()Provide default boundaries for optional arguments. fit_variogram
(x_data, y_data[, maxfev])Fiting the isotropic variogram-model to given data. fix_dim
()Set a fix dimension for the model. ln_spectral_rad_pdf
(r)Log radial spectral density of the model. percentile_scale
([per])Calculate the percentile scale of the isotrope model. plot
([func])Plot a function of a the CovModel. pykrige_vario
([args, r])Isotropic variogram of the model for pykrige. set_arg_bounds
(**kwargs)Set bounds for the parameters of the model. spectral_density
(k)Spectral density of the covariance model. spectral_rad_pdf
(r)Radial spectral density of the model. spectrum
(k)Spectrum of the covariance model. var_factor
()Factor for the variance. vario_nugget
(r)Isotropic variogram of the model respecting the nugget at r=0. vario_spatial
(pos)Spatial variogram respecting anisotropy and rotation. -
check_opt_arg
()[source]¶ Run checks for the optional arguments.
This is in addition to the bound-checks
Notes
- You can use this to raise a ValueError/warning
- Any return value will be ignored
- This method will only be run once, when the class is initialized
-
cov_nugget
(r)[source]¶ Covariance of the model respecting the nugget at r=0.
Given by:
Where is the correlation function.
-
default_opt_arg
()[source]¶ Provide default optional arguments by the user.
Should be given as a dictionary.
-
fit_variogram
(x_data, y_data, maxfev=1000, **para_deselect)[source]¶ Fiting the isotropic variogram-model to given data.
Parameters: - x_data (
numpy.ndarray
) – The radii of the meassured variogram. - y_data (
numpy.ndarray
) – The messured variogram - maxfev (int, optional) – The maximum number of calls to the function in scipy curvefit. Default: 1000
- **para_deselect – You can deselect the parameters to be fitted, by setting them “False” as keywords. By default, all parameters are fitted.
Returns: - fit_para (
dict
) – Dictonary with the fitted parameter values - pcov (
numpy.ndarray
) – The estimated covariance of popt fromscipy.optimize.curve_fit
Notes
You can set the bounds for each parameter by accessing
CovModel.set_arg_bounds
.The fitted parameters will be instantly set in the model.
- x_data (
-
percentile_scale
(per=0.9)[source]¶ Calculate the percentile scale of the isotrope model.
This is the distance, where the given percentile of the variance is reached by the variogram
-
plot
(func='variogram', **kwargs)[source]¶ Plot a function of a the CovModel.
Parameters: - func (
str
, optional) –Function to be plotted. Could be one of:
- ”variogram”
- ”covariance”
- ”correlation”
- ”vario_spatial”
- ”cov_spatial”
- ”cor_spatial”
- ”spectrum”
- ”spectral_density”
- ”spectral_rad_pdf”
- **kwargs – Keyword arguments forwarded to the plotting function
“plot_” + func in
gstools.covmodel.plot
.
See also
- func (
-
pykrige_vario
(args=None, r=0)[source]¶ Isotropic variogram of the model for pykrige.
Given by:
Where is the correlation function.
-
set_arg_bounds
(**kwargs)[source]¶ Set bounds for the parameters of the model.
Parameters: **kwargs – Parameter name as keyword (“var”, “len_scale”, “nugget”, <opt_arg>) and a list of 2 or 3 values as value:
[a, b]
or[a, b, <type>]
<type> is one of
"oo"
,"cc"
,"oc"
or"co"
to define if the bounds are open (“o”) or closed (“c”).
-
spectral_density
(k)[source]¶ Spectral density of the covariance model.
This is given by:
Where is the spectrum of the covariance model.
Parameters: k ( float
) – Radius of the phase:
-
spectrum
(k)[source]¶ Spectrum of the covariance model.
This is given by:
Internally, this is calculated by the hankel transformation:
Where is the covariance function of the model.
Parameters: k ( float
) – Radius of the phase:
-
vario_nugget
(r)[source]¶ Isotropic variogram of the model respecting the nugget at r=0.
Given by:
Where is the correlation function.
-
angles
¶ Rotation angles (in rad) of the model.
Type: numpy.ndarray
-
anis
¶ The anisotropy factors of the model.
Type: numpy.ndarray
-
arg_bounds
¶ Bounds for all parameters.
Notes
Keys are the opt-arg names and values are lists of 2 or 3 values:
[a, b]
or[a, b, <type>]
<type> is one of
"oo"
,"cc"
,"oc"
or"co"
to define if the bounds are open (“o”) or closed (“c”).Type: dict
-
hankel_kw
¶ hankel.SymmetricFourierTransform
kwargs.Type: dict
-
integral_scale
¶ The main integral scale of the model.
Raises: ValueError
– If integral scale is not setable.Type: float
-
integral_scale_vec
¶ The integral scales in each direction.
Notes
- This is calculated by:
integral_scale_vec[0] = integral_scale
integral_scale_vec[1] = integral_scale*anis[0]
integral_scale_vec[2] = integral_scale*anis[1]
Type: numpy.ndarray
-
len_scale_bounds
¶ Bounds for the lenght scale.
Notes
Is a list of 2 or 3 values:
[a, b]
or[a, b, <type>]
<type> is one of
"oo"
,"cc"
,"oc"
or"co"
to define if the bounds are open (“o”) or closed (“c”).Type: list
-
len_scale_vec
¶ The length scales in each direction.
Notes
- This is calculated by:
len_scale_vec[0] = len_scale
len_scale_vec[1] = len_scale*anis[0]
len_scale_vec[2] = len_scale*anis[1]
Type: numpy.ndarray
-
nugget_bounds
¶ Bounds for the nugget.
Notes
Is a list of 2 or 3 values:
[a, b]
or[a, b, <type>]
<type> is one of
"oo"
,"cc"
,"oc"
or"co"
to define if the bounds are open (“o”) or closed (“c”).Type: list
-
opt_arg_bounds
¶ Bounds for the optional arguments.
Notes
Keys are the opt-arg names and values are lists of 2 or 3 values:
[a, b]
or[a, b, <type>]
<type> is one of
"oo"
,"cc"
,"oc"
or"co"
to define if the bounds are open (“o”) or closed (“c”).Type: dict
-
pykrige_angle
¶ 2D rotation angle for pykrige.
-
pykrige_angle_x
¶ 3D rotation angle around x for pykrige.
-
pykrige_angle_y
¶ 3D rotation angle around y for pykrige.
-
pykrige_angle_z
¶ 3D rotation angle around z for pykrige.
-
pykrige_anis
¶ 2D anisotropy ratio for pykrige.
-
pykrige_anis_y
¶ 3D anisotropy ratio in y direction for pykrige.
-
pykrige_anis_z
¶ 3D anisotropy ratio in z direction for pykrige.
-
pykrige_kwargs
¶ Keyword arguments for pykrige routines.
- dim (