gstools.covmodel.models

GStools subpackage providing different covariance models.

The following classes and functions are provided

Gaussian([dim, var, len_scale, nugget, …]) The Gaussian covariance model
Exponential([dim, var, len_scale, nugget, …]) The Exponential covariance model
Matern([dim, var, len_scale, nugget, anis, …]) The Matérn covariance model
Rational([dim, var, len_scale, nugget, …]) The rational quadratic covariance model
Stable([dim, var, len_scale, nugget, anis, …]) The stable covariance model
Spherical([dim, var, len_scale, nugget, …]) The Spherical covariance model
Linear([dim, var, len_scale, nugget, anis, …]) The bounded linear covariance model
MaternRescal([dim, var, len_scale, nugget, …]) The rescaled Matérn covariance model
SphericalRescal([dim, var, len_scale, …]) The rescaled Spherical covariance model
class gstools.covmodel.models.Gaussian(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: gstools.covmodel.base.CovModel

The Gaussian covariance model

Notes

This model is given by the following correlation function:

\mathrm{cor}(r) =
\exp\left(- \frac{\pi}{4} \cdot \left(\frac{r}{\ell}\right)^2\right)

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 or list, 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 or list, optional) – anisotropy ratios in the transversal directions [y, z]. Default: 1.0
  • angles (float or list, 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 or list or None, 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 or None, optional) – raw variance of the model which will be multiplied with CovModel.var_factor to result in the actual variance. If given, var will be ignored. (This is just for models that override CovModel.var_factor) Default: None
  • hankel_kw (dict or None, optional) – Modify the init-arguments of hankel.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

Methods

calc_integral_scale() The integral scale of the gaussian model is the length scale
correlation(r) Gaussian correlation function
covariance(r) Covariance of the model
spectral_rad_cdf(r) The cdf of the radial spectral density
spectral_rad_ppf(u) The ppf of the radial spectral density
spectrum(k) The spectrum of the covariance model.
variogram(r) Isotropic variogram of the model
variogram_normed(r) Normalized variogram of the model
calc_integral_scale()[source]

The integral scale of the gaussian model is the length scale

correlation(r)[source]

Gaussian correlation function

\mathrm{cor}(r) =
\exp\left(- \frac{\pi}{4} \cdot \left(\frac{r}{\ell}\right)^2\right)

covariance(r)

Covariance of the model

Given by: C\left(r\right)=
\sigma^2\cdot\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

spectral_rad_cdf(r)[source]

The cdf of the radial spectral density

spectral_rad_ppf(u)[source]

The ppf of the radial spectral density

spectrum(k)[source]

The spectrum of the covariance model.

This is given by:

S(k) = \left(\frac{1}{2\pi}\right)^n
\int C(r) e^{i b\mathbf{k}\cdot\mathbf{r}} d^n\mathbf{r}

Internally, this is calculated by the hankel transformation:

S(k) = \left(\frac{1}{2\pi}\right)^n \cdot
\frac{(2\pi)^{n/2}}{(bk)^{n/2-1}}
\int_0^\infty r^{n/2-1} f(r) J_{n/2-1}(bkr) r dr

Where C(r) is the covariance function of the model.

Parameters:k (float) – Radius of the phase: k=\left\Vert\mathbf{k}\right\Vert
variogram(r)

Isotropic variogram of the model

Given by: \gamma\left(r\right)=
\sigma^2\cdot\left(1-\mathrm{cor}\left(r\right)\right)+n

Where \mathrm{cor}(r) is the correlation function.

variogram_normed(r)

Normalized variogram of the model

Given by: \tilde{\gamma}\left(r\right)=
1-\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

class gstools.covmodel.models.Exponential(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: gstools.covmodel.base.CovModel

The Exponential covariance model

Notes

This model is given by the following correlation function:

\mathrm{cor}(r) =
\exp\left(- \frac{r}{\ell} \right)

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 or list, 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 or list, optional) – anisotropy ratios in the transversal directions [y, z]. Default: 1.0
  • angles (float or list, 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 or list or None, 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 or None, optional) – raw variance of the model which will be multiplied with CovModel.var_factor to result in the actual variance. If given, var will be ignored. (This is just for models that override CovModel.var_factor) Default: None
  • hankel_kw (dict or None, optional) – Modify the init-arguments of hankel.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

Methods

calc_integral_scale() The integral scale of the exponential model is the length scale
correlation(r) Exponential correlation function
covariance(r) Covariance of the model
spectral_rad_cdf(r) The cdf of the radial spectral density
spectral_rad_ppf(u) The ppf of the radial spectral density
spectrum(k) The spectrum of the covariance model.
variogram(r) Isotropic variogram of the model
variogram_normed(r) Normalized variogram of the model
calc_integral_scale()[source]

The integral scale of the exponential model is the length scale

correlation(r)[source]

Exponential correlation function

\mathrm{cor}(r) =
\exp\left(- \frac{r}{\ell} \right)

covariance(r)

Covariance of the model

Given by: C\left(r\right)=
\sigma^2\cdot\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

spectral_rad_cdf(r)[source]

The cdf of the radial spectral density

spectral_rad_ppf(u)[source]

The ppf of the radial spectral density

spectrum(k)[source]

The spectrum of the covariance model.

This is given by:

S(k) = \left(\frac{1}{2\pi}\right)^n
\int C(r) e^{i b\mathbf{k}\cdot\mathbf{r}} d^n\mathbf{r}

Internally, this is calculated by the hankel transformation:

S(k) = \left(\frac{1}{2\pi}\right)^n \cdot
\frac{(2\pi)^{n/2}}{(bk)^{n/2-1}}
\int_0^\infty r^{n/2-1} f(r) J_{n/2-1}(bkr) r dr

Where C(r) is the covariance function of the model.

Parameters:k (float) – Radius of the phase: k=\left\Vert\mathbf{k}\right\Vert
variogram(r)

Isotropic variogram of the model

Given by: \gamma\left(r\right)=
\sigma^2\cdot\left(1-\mathrm{cor}\left(r\right)\right)+n

Where \mathrm{cor}(r) is the correlation function.

variogram_normed(r)

Normalized variogram of the model

Given by: \tilde{\gamma}\left(r\right)=
1-\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

class gstools.covmodel.models.Spherical(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: gstools.covmodel.base.CovModel

The Spherical covariance model

Notes

This model is given by the following correlation function:

\mathrm{cor}(r) =
\begin{cases}
1-\frac{3}{2}\cdot\frac{r}{\ell} +
\frac{1}{2}\cdot\left(\frac{r}{\ell}\right)^{3}
& r<\ell\\
0 & r\geq\ell
\end{cases}

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 or list, 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 or list, optional) – anisotropy ratios in the transversal directions [y, z]. Default: 1.0
  • angles (float or list, 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 or list or None, 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 or None, optional) – raw variance of the model which will be multiplied with CovModel.var_factor to result in the actual variance. If given, var will be ignored. (This is just for models that override CovModel.var_factor) Default: None
  • hankel_kw (dict or None, optional) – Modify the init-arguments of hankel.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

Methods

correlation(r) Spherical correlation function
covariance(r) Covariance of the model
variogram(r) Isotropic variogram of the model
variogram_normed(r) Normalized variogram of the model
correlation(r)[source]

Spherical correlation function

\mathrm{cor}(r) =
\begin{cases}
1-\frac{3}{2}\cdot\frac{r}{\ell} +
\frac{1}{2}\cdot\left(\frac{r}{\ell}\right)^{3}
& r<\ell\\
0 & r\geq\ell
\end{cases}

covariance(r)

Covariance of the model

Given by: C\left(r\right)=
\sigma^2\cdot\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

variogram(r)

Isotropic variogram of the model

Given by: \gamma\left(r\right)=
\sigma^2\cdot\left(1-\mathrm{cor}\left(r\right)\right)+n

Where \mathrm{cor}(r) is the correlation function.

variogram_normed(r)

Normalized variogram of the model

Given by: \tilde{\gamma}\left(r\right)=
1-\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

class gstools.covmodel.models.SphericalRescal(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: gstools.covmodel.base.CovModel

The rescaled Spherical covariance model

Notes

This model is given by the following correlation function:

\mathrm{cor}(r) =
\begin{cases}
1-\frac{9}{16}\cdot\frac{r}{\ell} +
\frac{27}{1024}\cdot\left(\frac{r}{\ell}\right)^{3}
& r<\frac{8}{3}\ell\\
0 & r\geq\frac{8}{3}\ell
\end{cases}

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 or list, 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 or list, optional) – anisotropy ratios in the transversal directions [y, z]. Default: 1.0
  • angles (float or list, 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 or list or None, 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 or None, optional) – raw variance of the model which will be multiplied with CovModel.var_factor to result in the actual variance. If given, var will be ignored. (This is just for models that override CovModel.var_factor) Default: None
  • hankel_kw (dict or None, optional) – Modify the init-arguments of hankel.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

Methods

calc_integral_scale() The integral scale of this spherical model is the length scale
correlation(r) Rescaled Spherical correlation function
covariance(r) Covariance of the model
variogram(r) Isotropic variogram of the model
variogram_normed(r) Normalized variogram of the model
calc_integral_scale()[source]

The integral scale of this spherical model is the length scale

correlation(r)[source]

Rescaled Spherical correlation function

\mathrm{cor}(r) =
\begin{cases}
1-\frac{9}{16}\cdot\frac{r}{\ell} +
\frac{27}{1024}\cdot\left(\frac{r}{\ell}\right)^{3}
& r<\frac{8}{3}\ell\\
0 & r\geq\frac{8}{3}\ell
\end{cases}

covariance(r)

Covariance of the model

Given by: C\left(r\right)=
\sigma^2\cdot\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

variogram(r)

Isotropic variogram of the model

Given by: \gamma\left(r\right)=
\sigma^2\cdot\left(1-\mathrm{cor}\left(r\right)\right)+n

Where \mathrm{cor}(r) is the correlation function.

variogram_normed(r)

Normalized variogram of the model

Given by: \tilde{\gamma}\left(r\right)=
1-\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

class gstools.covmodel.models.Rational(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: gstools.covmodel.base.CovModel

The rational quadratic covariance model

Notes

This model is given by the following correlation function:

\mathrm{cor}(r) =
\left(1 + \frac{1}{2\alpha} \cdot
\left(\frac{r}{\ell}\right)^2\right)^{-\alpha}

\alpha is a shape parameter and should be > 0.5.

Other Parameters:
 
  • **opt_arg – The following parameters are covered by these keyword arguments
  • alpha (float, optional) – Shape parameter. Standard range: (0, inf) Default: 1.0
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 or list, 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 or list, optional) – anisotropy ratios in the transversal directions [y, z]. Default: 1.0
  • angles (float or list, 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 or list or None, 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 or None, optional) – raw variance of the model which will be multiplied with CovModel.var_factor to result in the actual variance. If given, var will be ignored. (This is just for models that override CovModel.var_factor) Default: None
  • hankel_kw (dict or None, optional) – Modify the init-arguments of hankel.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

Methods

correlation(r) Rational correlation function
covariance(r) Covariance of the model
default_opt_arg() The defaults for the optional arguments:
default_opt_arg_bounds() The defaults boundaries for the optional arguments:
variogram(r) Isotropic variogram of the model
variogram_normed(r) Normalized variogram of the model
correlation(r)[source]

Rational correlation function

\mathrm{cor}(r) =
\left(1 + \frac{1}{2\alpha} \cdot
\left(\frac{r}{\ell}\right)^2\right)^{-\alpha}

covariance(r)

Covariance of the model

Given by: C\left(r\right)=
\sigma^2\cdot\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

default_opt_arg()[source]

The defaults for the optional arguments:

  • {"alpha": 1.0}
Returns:Defaults for optional arguments
Return type:dict
default_opt_arg_bounds()[source]

The defaults boundaries for the optional arguments:

  • {"alpha": [0.5, inf]}
Returns:Boundaries for optional arguments
Return type:dict
variogram(r)

Isotropic variogram of the model

Given by: \gamma\left(r\right)=
\sigma^2\cdot\left(1-\mathrm{cor}\left(r\right)\right)+n

Where \mathrm{cor}(r) is the correlation function.

variogram_normed(r)

Normalized variogram of the model

Given by: \tilde{\gamma}\left(r\right)=
1-\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

class gstools.covmodel.models.Stable(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: gstools.covmodel.base.CovModel

The stable covariance model

Notes

This model is given by the following correlation function:

\mathrm{cor}(r) =
\exp\left(- \left(\frac{r}{\ell}\right)^{\alpha}\right)

\alpha is a shape parameter with \alpha\in(0,2]

Other Parameters:
 
  • **opt_arg – The following parameters are covered by these keyword arguments
  • alpha (float, optional) – Shape parameter. Standard range: (0, 2] Default: 1.5
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 or list, 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 or list, optional) – anisotropy ratios in the transversal directions [y, z]. Default: 1.0
  • angles (float or list, 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 or list or None, 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 or None, optional) – raw variance of the model which will be multiplied with CovModel.var_factor to result in the actual variance. If given, var will be ignored. (This is just for models that override CovModel.var_factor) Default: None
  • hankel_kw (dict or None, optional) – Modify the init-arguments of hankel.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

Methods

check_opt_arg() Checks for the optional arguments
correlation(r) Stable correlation function
covariance(r) Covariance of the model
default_opt_arg() The defaults for the optional arguments:
default_opt_arg_bounds() The defaults boundaries for the optional arguments:
variogram(r) Isotropic variogram of the model
variogram_normed(r) Normalized variogram of the model
check_opt_arg()[source]

Checks for the optional arguments

Warns:alpha – If alpha is < 0.3, the model tends to a nugget model and gets numerically unstable.
correlation(r)[source]

Stable correlation function

\mathrm{cor}(r) =
\exp\left(- \left(\frac{r}{\ell}\right)^{\alpha}\right)

covariance(r)

Covariance of the model

Given by: C\left(r\right)=
\sigma^2\cdot\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

default_opt_arg()[source]

The defaults for the optional arguments:

  • {"alpha": 1.5}
Returns:Defaults for optional arguments
Return type:dict
default_opt_arg_bounds()[source]

The defaults boundaries for the optional arguments:

  • {"alpha": [0, 2, "oc"]}
Returns:Boundaries for optional arguments
Return type:dict
variogram(r)

Isotropic variogram of the model

Given by: \gamma\left(r\right)=
\sigma^2\cdot\left(1-\mathrm{cor}\left(r\right)\right)+n

Where \mathrm{cor}(r) is the correlation function.

variogram_normed(r)

Normalized variogram of the model

Given by: \tilde{\gamma}\left(r\right)=
1-\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

class gstools.covmodel.models.Matern(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: gstools.covmodel.base.CovModel

The Matérn covariance model

Notes

This model is given by the following correlation function:

\mathrm{cor}(r) =
\frac{2^{1-\nu}}{\Gamma\left(\nu\right)} \cdot
\left(\sqrt{2\nu}\cdot\frac{r}{\ell}\right)^{\nu} \cdot
\mathrm{K}_{\nu}\left(\sqrt{2\nu}\cdot\frac{r}{\ell}\right)

Where \Gamma is the gamma function and \mathrm{K}_{\nu} is the modified Bessel function of the second kind.

\nu is a shape parameter and should be >= 0.5.

Other Parameters:
 
  • **opt_arg – The following parameters are covered by these keyword arguments
  • nu (float, optional) – Shape parameter. Standard range: [0.5, 60] Default: 1.0
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 or list, 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 or list, optional) – anisotropy ratios in the transversal directions [y, z]. Default: 1.0
  • angles (float or list, 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 or list or None, 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 or None, optional) – raw variance of the model which will be multiplied with CovModel.var_factor to result in the actual variance. If given, var will be ignored. (This is just for models that override CovModel.var_factor) Default: None
  • hankel_kw (dict or None, optional) – Modify the init-arguments of hankel.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

Methods

check_opt_arg() Checks for the optional arguments
correlation(r) Matérn correlation function
covariance(r) Covariance of the model
default_opt_arg() The defaults for the optional arguments:
default_opt_arg_bounds() The defaults boundaries for the optional arguments:
variogram(r) Isotropic variogram of the model
variogram_normed(r) Normalized variogram of the model
check_opt_arg()[source]

Checks for the optional arguments

Warns:nu – If nu is > 50, the model gets numerically unstable.
correlation(r)[source]

Matérn correlation function

\mathrm{cor}(r) =
\frac{2^{1-\nu}}{\Gamma\left(\nu\right)} \cdot
\left(\sqrt{2\nu}\cdot\frac{r}{\ell}\right)^{\nu} \cdot
\mathrm{K}_{\nu}\left(\sqrt{2\nu}\cdot\frac{r}{\ell}\right)

covariance(r)

Covariance of the model

Given by: C\left(r\right)=
\sigma^2\cdot\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

default_opt_arg()[source]

The defaults for the optional arguments:

  • {"nu": 1.0}
Returns:Defaults for optional arguments
Return type:dict
default_opt_arg_bounds()[source]

The defaults boundaries for the optional arguments:

  • {"nu": [0.5, 60.0, "cc"]}
Returns:Boundaries for optional arguments
Return type:dict
variogram(r)

Isotropic variogram of the model

Given by: \gamma\left(r\right)=
\sigma^2\cdot\left(1-\mathrm{cor}\left(r\right)\right)+n

Where \mathrm{cor}(r) is the correlation function.

variogram_normed(r)

Normalized variogram of the model

Given by: \tilde{\gamma}\left(r\right)=
1-\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

class gstools.covmodel.models.MaternRescal(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: gstools.covmodel.base.CovModel

The rescaled Matérn covariance model

Notes

This model is given by the following correlation function:

\mathrm{cor}(r) =
\frac{2^{1-\nu}}{\Gamma\left(\nu\right)} \cdot
\left(\frac{\pi}{B\left(\nu,\frac{1}{2}\right)} \cdot
\frac{r}{\ell}\right)^{\nu} \cdot
\mathrm{K}_{\nu}\left(\frac{\pi}{B\left(\nu,\frac{1}{2}\right)} \cdot
\frac{r}{\ell}\right)

Where \Gamma is the gamma function, \mathrm{K}_{\nu} is the modified Bessel function of the second kind and B is the Euler beta function.

\nu is a shape parameter and should be > 0.5.

Other Parameters:
 
  • **opt_arg – The following parameters are covered by these keyword arguments
  • nu (float, optional) – Shape parameter. Standard range: [0.5, 60] Default: 1.0
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 or list, 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 or list, optional) – anisotropy ratios in the transversal directions [y, z]. Default: 1.0
  • angles (float or list, 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 or list or None, 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 or None, optional) – raw variance of the model which will be multiplied with CovModel.var_factor to result in the actual variance. If given, var will be ignored. (This is just for models that override CovModel.var_factor) Default: None
  • hankel_kw (dict or None, optional) – Modify the init-arguments of hankel.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

Methods

check_opt_arg() Checks for the optional arguments
correlation(r) Rescaled Matérn correlation function
covariance(r) Covariance of the model
default_opt_arg() The defaults for the optional arguments:
default_opt_arg_bounds() The defaults boundaries for the optional arguments:
variogram(r) Isotropic variogram of the model
variogram_normed(r) Normalized variogram of the model
check_opt_arg()[source]

Checks for the optional arguments

Warns:nu – If nu is > 50, the model gets numerically unstable.
correlation(r)[source]

Rescaled Matérn correlation function

\mathrm{cor}(r) =
\frac{2^{1-\nu}}{\Gamma\left(\nu\right)} \cdot
\left(\frac{\pi}{B\left(\nu,\frac{1}{2}\right)} \cdot
\frac{r}{\ell}\right)^{\nu} \cdot
\mathrm{K}_{\nu}\left(\frac{\pi}{B\left(\nu,\frac{1}{2}\right)}
\cdot\frac{r}{\ell}\right)

covariance(r)

Covariance of the model

Given by: C\left(r\right)=
\sigma^2\cdot\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

default_opt_arg()[source]

The defaults for the optional arguments:

  • {"nu": 1.0}
Returns:Defaults for optional arguments
Return type:dict
default_opt_arg_bounds()[source]

The defaults boundaries for the optional arguments:

  • {"nu": [0.5, 60.0, "cc"]}
Returns:Boundaries for optional arguments
Return type:dict
variogram(r)

Isotropic variogram of the model

Given by: \gamma\left(r\right)=
\sigma^2\cdot\left(1-\mathrm{cor}\left(r\right)\right)+n

Where \mathrm{cor}(r) is the correlation function.

variogram_normed(r)

Normalized variogram of the model

Given by: \tilde{\gamma}\left(r\right)=
1-\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

class gstools.covmodel.models.Linear(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: gstools.covmodel.base.CovModel

The bounded linear covariance model

Notes

This model is given by the following correlation function:

\mathrm{cor}(r) =
\begin{cases}
1-\frac{r}{\ell}
& r<\ell\\
0 & r\geq\ell
\end{cases}

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 or list, 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 or list, optional) – anisotropy ratios in the transversal directions [y, z]. Default: 1.0
  • angles (float or list, 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 or list or None, 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 or None, optional) – raw variance of the model which will be multiplied with CovModel.var_factor to result in the actual variance. If given, var will be ignored. (This is just for models that override CovModel.var_factor) Default: None
  • hankel_kw (dict or None, optional) – Modify the init-arguments of hankel.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

Methods

correlation(r) Linear correlation function
covariance(r) Covariance of the model
variogram(r) Isotropic variogram of the model
variogram_normed(r) Normalized variogram of the model
correlation(r)[source]

Linear correlation function

\mathrm{cor}(r) =
\begin{cases}
1-\frac{r}{\ell}
& r<\ell\\
0 & r\geq\ell
\end{cases}

covariance(r)

Covariance of the model

Given by: C\left(r\right)=
\sigma^2\cdot\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.

variogram(r)

Isotropic variogram of the model

Given by: \gamma\left(r\right)=
\sigma^2\cdot\left(1-\mathrm{cor}\left(r\right)\right)+n

Where \mathrm{cor}(r) is the correlation function.

variogram_normed(r)

Normalized variogram of the model

Given by: \tilde{\gamma}\left(r\right)=
1-\mathrm{cor}\left(r\right)

Where \mathrm{cor}(r) is the correlation function.