gstools.covmodel.tpl_models¶
GStools subpackage providing truncated power law covariance models.
The following classes and functions are provided
TPLGaussian ([dim, var, len_scale, nugget, …]) |
Truncated-Power-Law with Gaussian modes. |
TPLExponential ([dim, var, len_scale, …]) |
Truncated-Power-Law with Exponential modes. |
TPLStable ([dim, var, len_scale, nugget, …]) |
Truncated-Power-Law with Stable modes. |
-
class
gstools.covmodel.tpl_models.
TPLGaussian
(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
Truncated-Power-Law with Gaussian modes.
Notes
The truncated power law is given by a superposition of scale-dependent variograms:
with Gaussian modes on each scale:
This results in:
The “length scale” of this model is equivalent by the integration range:
If you want to define an upper scale truncation, you should set
len_low
andlen_scale
accordingly.The following Parameters occure:
: The scaling factor from the Power-Law. This parameter will be calculated internally by the given variance. You can access C directly by
model.var_raw
: The hurst coefficient (
model.hurst
): The lower length scale truncation of the model (
model.len_low
): The upper length scale truncation of the model (
model.len_up
)This will be calculated internally by:
len_up = len_low + len_scale
That means, that the
len_scale
in this model actually represents the integration range for the truncated power law.is the exponential integral.
Other Parameters: 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}
.
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.len_up
float
: Upper length scale truncation of the model.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. correlation
(r)Truncated-Power-Law with Gaussian modes - correlation function. cov_nugget
(r)Covariance of the model respecting the nugget at r=0. cov_spatial
(pos)Spatial covariance respecting anisotropy and rotation. covariance
(r)Covariance of the model. default_arg_bounds
()Provide default boundaries for arguments. default_opt_arg
()Defaults for the optional arguments. default_opt_arg_bounds
()Defaults for boundaries of the 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 C (Power-Law factor) to result in variance. vario_nugget
(r)Isotropic variogram of the model respecting the nugget at r=0. vario_spatial
(pos)Spatial variogram respecting anisotropy and rotation. variogram
(r)Isotropic variogram of the model. -
correlation
(r)[source]¶ Truncated-Power-Law with Gaussian modes - correlation function.
If
len_low=0
we have a simple representation:The general case:
-
covariance
(r)¶ Covariance of the model.
Given by:
Where is the correlation function.
-
default_opt_arg
()[source]¶ Defaults for the optional arguments.
{"hurst": 0.5, "len_low": 0.0}
Returns: Defaults for optional arguments Return type: dict
-
default_opt_arg_bounds
()[source]¶ Defaults for boundaries of the optional arguments.
{"hurst": [0, 1, "oo"], "len_low": [0, 1000, "cc"]}
Returns: Boundaries for optional arguments Return type: dict
-
var_factor
()[source]¶ Factor for C (Power-Law factor) to result in variance.
This is used to result in the right variance, which is depending on the hurst coefficient and the length-scale extents
Returns: factor Return type: float
-
variogram
(r)¶ Isotropic variogram of the model.
Given by:
Where is the correlation function.
-
class
gstools.covmodel.tpl_models.
TPLExponential
(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
Truncated-Power-Law with Exponential modes.
Notes
The truncated power law is given by a superposition of scale-dependent variograms:
with Exponential modes on each scale:
This results in:
The “length scale” of this model is equivalent by the integration range:
If you want to define an upper scale truncation, you should set
len_low
andlen_scale
accordingly.The following Parameters occure:
: The scaling factor from the Power-Law. This parameter will be calculated internally by the given variance. You can access C directly by
model.var_raw
: The hurst coefficient (
model.hurst
): The lower length scale truncation of the model (
model.len_low
): The upper length scale truncation of the model (
model.len_up
)This will be calculated internally by:
len_up = len_low + len_scale
That means, that the
len_scale
in this model actually represents the integration range for the truncated power law.is the exponential integral.
Other Parameters: 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}
.
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.len_up
float
: Upper length scale truncation of the model.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. correlation
(r)Truncated-Power-Law with Exponential modes - correlation function. cov_nugget
(r)Covariance of the model respecting the nugget at r=0. cov_spatial
(pos)Spatial covariance respecting anisotropy and rotation. covariance
(r)Covariance of the model. default_arg_bounds
()Provide default boundaries for arguments. default_opt_arg
()Defaults for the optional arguments. default_opt_arg_bounds
()Defaults for boundaries of the 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 C (Power-Law factor) to result in variance. vario_nugget
(r)Isotropic variogram of the model respecting the nugget at r=0. vario_spatial
(pos)Spatial variogram respecting anisotropy and rotation. variogram
(r)Isotropic variogram of the model. -
correlation
(r)[source]¶ Truncated-Power-Law with Exponential modes - correlation function.
If
len_low=0
we have a simple representation:The general case:
-
covariance
(r)¶ Covariance of the model.
Given by:
Where is the correlation function.
-
default_opt_arg
()[source]¶ Defaults for the optional arguments.
{"hurst": 0.25, "len_low": 0.0}
Returns: Defaults for optional arguments Return type: dict
-
default_opt_arg_bounds
()[source]¶ Defaults for boundaries of the optional arguments.
{"hurst": [0, 1, "oo"], "len_low": [0, 1000, "cc"]}
Returns: Boundaries for optional arguments Return type: dict
-
var_factor
()[source]¶ Factor for C (Power-Law factor) to result in variance.
This is used to result in the right variance, which is depending on the hurst coefficient and the length-scale extents
Returns: factor Return type: float
-
variogram
(r)¶ Isotropic variogram of the model.
Given by:
Where is the correlation function.
-
class
gstools.covmodel.tpl_models.
TPLStable
(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
Truncated-Power-Law with Stable modes.
Notes
The truncated power law is given by a superposition of scale-dependent variograms:
with Stable modes on each scale:
This results in:
The “length scale” of this model is equivalent by the integration range:
If you want to define an upper scale truncation, you should set
len_low
andlen_scale
accordingly.The following Parameters occure:
: The shape parameter of the Stable model.
- : Exponential modes
- : Gaussian modes
: The scaling factor from the Power-Law. This parameter will be calculated internally by the given variance. You can access C directly by
model.var_raw
: The hurst coefficient (
model.hurst
): The lower length scale truncation of the model (
model.len_low
): The upper length scale truncation of the model (
model.len_up
)This will be calculated internally by:
len_up = len_low + len_scale
That means, that the
len_scale
in this model actually represents the integration range for the truncated power law.is the exponential integral.
Other Parameters: - **opt_arg – The following parameters are covered by these keyword arguments
- hurst (
float
, optional) – Hurst coefficient of the power law. Standard range:(0, 1)
. Default:0.5
- alpha (
float
, optional) – Shape parameter of the stable model. Standard range:(0, 2]
. Default:1.5
- len_low (
float
, optional) – The lower length scale truncation of the model. Standard range:[0, 1000]
. Default:0.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
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}
.
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.len_up
float
: Upper length scale truncation of the model.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
()Check the optional arguments. cor_spatial
(pos)Spatial correlation respecting anisotropy and rotation. correlation
(r)Truncated-Power-Law with Stable modes - correlation function. cov_nugget
(r)Covariance of the model respecting the nugget at r=0. cov_spatial
(pos)Spatial covariance respecting anisotropy and rotation. covariance
(r)Covariance of the model. default_arg_bounds
()Provide default boundaries for arguments. default_opt_arg
()Defaults for the optional arguments. default_opt_arg_bounds
()Defaults for boundaries of the 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 C (Power-Law factor) to result in variance. vario_nugget
(r)Isotropic variogram of the model respecting the nugget at r=0. vario_spatial
(pos)Spatial variogram respecting anisotropy and rotation. variogram
(r)Isotropic variogram of the model. -
check_opt_arg
()[source]¶ Check the optional arguments.
Warns: alpha – If alpha is < 0.3, the model tends to a nugget model and gets numerically unstable.
-
correlation
(r)[source]¶ Truncated-Power-Law with Stable modes - correlation function.
If
len_low=0
we have a simple representation:The general case:
-
covariance
(r)¶ Covariance of the model.
Given by:
Where is the correlation function.
-
default_opt_arg
()[source]¶ Defaults for the optional arguments.
{"hurst": 0.5, "alpha": 1.5, "len_low": 0.0}
Returns: Defaults for optional arguments Return type: dict
-
default_opt_arg_bounds
()[source]¶ Defaults for boundaries of the optional arguments.
{"hurst": [0, 1, "oo"], "alpha": [0, 2, "oc"], "len_low": [0, 1000, "cc"]}
Returns: Boundaries for optional arguments Return type: dict
-
var_factor
()[source]¶ Factor for C (Power-Law factor) to result in variance.
This is used to result in the right variance, which is depending on the hurst coefficient and the length-scale extents
Returns: factor Return type: float
-
variogram
(r)¶ Isotropic variogram of the model.
Given by:
Where is the correlation function.