spam.excursions package#
Submodules#
spam.excursions.elkc module#
Library of SPAM functions for computation Lipschitz-Killing Curvatures. Copyright (C) 2020 SPAM Contributors
This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.
This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.
- spam.excursions.elkc.expectedMesures(kappa, j, n, hittingSet='tail', distributionType='gauss', mu=0.0, std=1.0, correlationType='gauss', lc=1.0, nu=2.0, a=1.0)[source]#
Compute the Lipschitz-Killing Curvatures E(LKC)
- Parameters:
kappa (float or list) – value of the threshold
j (int) – number of the functionnal
n (int) – spatial dimension
hittingSet (string) – the hitting set
distributionType (string) – type of distribution (see gaussianMinkowskiFunctionals)
mu (float) – mean value of the distribution
std (float) – standard deviation of the distribution
correlationType (string) – type of correlation function
lc (float) – correlation length
nu (float) – parameter used for correlationType=’matern’
a (float or list of floats) – size(s) of the object if float, the object is considered as a cube
- Returns:
the Gaussian Minkowski functionnal
- Return type:
Same type as kappa
- spam.excursions.elkc.gaussianMinkowskiFunctionals(kappa, j, hittingSet='tail', distributionType='gauss', mu=0.0, std=1.0)[source]#
Evaluate the Gaussian Minkowski functionals j
- Parameters:
kappa (float or list of floats) – value(s) of the threshold(s)
j (int) – number of the functionnal
hittingSet (string) – hitting set
distributionType (string) – type of the underlying distribution. Implemented:
gauss
andlog
mu (float) – mean value of the distribution
std (float) – standard deviation of the distribution
- Returns:
mink – the Gaussian Minkowski functionnal
- Return type:
same type as kappa
- spam.excursions.elkc.secondSpectralMoment(std, lc, correlationType='gauss', nu=2.0)[source]#
Compute the second spectral moment for a given covariance function:
\[{\lambda}_2 = C''(h)|_{h=0}\]- Parameters:
std (float) – standard deviation
lc (float) – correlation length
correlationType (string) – type of correlation function
nu (float) – parameter used for correlationType=’matern’
- Returns:
the second spectral moment
- Return type:
float
- spam.excursions.elkc.flag(n, j)[source]#
Compute the flag coefficients
[ n, j ] = ( n, j ) w_n / w_{n-j}w_j
- Parameters:
n (int) – first parameter of the binom
j (int) – second parameter of the binom
- Returns:
the flag coefficient [n , j]
- Return type:
float
- spam.excursions.elkc.hermitePolynomial(x, n)[source]#
Evaluate a x the probabilitic Hermite polynomial n
- Parameters:
x (float) – point of evaluation
n (int) – number of Hermite polynomia
- Returns:
He_n(x)
the value of the polynomial- Return type:
float
- spam.excursions.elkc.ballVolume(n, r=1.0)[source]#
Compute the volume of a nth dimensional ball
- Parameters:
n (int) – spatial dimension
r (float) – radius of the ball
- Returns:
volume of the ball
- Return type:
float
- spam.excursions.elkc.lkcCuboid(i, n, a)[source]#
Compute the Lipschitz-Killing Curvatures of a parallelepiped
- Parameters:
i (int) – number of the LKC.
0 <= i <= n
n (int) – spatial dimension
a (float or list of float) – size(s) of the cuboid. If float, the object is considered to be a cube.
- Returns:
The Lipschitz-Killing Curvature
- Return type:
float
spam.excursions.randomFields module#
This module is a wrapper around the python library gstools made to intuitively generate realisations of correlated Random Fields on a regular grid with a given size, discretisation and characteristic length.
- spam.excursions.randomFields.simulateRandomField(lengths=1.0, nNodes=100, covarianceModel='Gaussian', covarianceParameters={}, dim=3, nRea=1, seed=None, vtkFile=None)[source]#
Generates realisations of Correlated Random Field based on the python library gstools.
- Parameters:
lengths (float | list of floats, default=1.0) – Length of the definition domain in every directions (should be size of dim). If a single float is given it will be used in each directions.
nNodes (int | list of int, default=100) – Number of nodes in every directions. If a single float is given it will be used in each directions.
covarianceModel (string, default="Gaussian") – The name of the covariance model used by gstools. To be chosen among this list of models.
covarianceParameters (dict, default={}) – A dictionnary of the keyword arguments of the covariance model used (see gstools documentation).
dim (int default=3) – Spatial dimention
nRea (int, default=1) – number of realisations
seed (int, default=None) – Defines the seed for the random numbers generator. Setting a seed to a certain integer yields the same realisations.
vtkFile (string, default=None) – If not None, the base name of the vtk files saved.
- Returns:
The realisations of the random fields.
- Return type:
List of nRea numpy arrays if shape nNodes
Example
>>> from spam.excursions import simulateRandomField >>> >>> # generation of two realisations of a Gaussian Correlated Random field >>> realisations = simulateRandomField( >>> lengths=25, >>> nNodes=100, >>> dim=3, >>> nRea=2, >>> covarianceModel="Gaussian", >>> covarianceParameters={"len_scale": 10}, >>> vtkFile="gaussian_randomfield" >>> ) >>> >>> for realisation in realisations: >>> print(realisation.shape) >>> >>> # generation of one realisation of a Matern Correlated Random field >>> realisations = simulateRandomField( >>> lengths=25, >>> nNodes=100, >>> dim=3, >>> covarianceModel="Matern", >>> covarianceParameters={"len_scale": 10, "nu": 0.2}, >>> vtkFile="matern_randomfield" >>> )
Warning
If no rescaling factors are given in the covariance parameters it is set to 1 which deviates from gstools default behavior.
- spam.excursions.randomFields.parametersLogToGauss(muLog, stdLog, lcLog=1)[source]#
Gives the underlying Gaussian standard deviation and mean value of the log normal distribution of standard deviation and mean value.
- Parameters:
muLog (float) – Mean value of the log normal distribution.
stdLog (float) – Standard deviation of the log normal distribution.
lcLog (float, default=1) – Correlation length of the underlying log normal correlated Random Field.
- Returns:
muGauss (float) – Mean value of the Gaussian distribution.
stdGauss (float) – Standard deviation of the Gaussian distribution.
lcGauss (float) – Correlation length of the underlying Gaussian correlated Random Field.
- spam.excursions.randomFields.fieldGaussToUniform(g, a=0.0, b=1.0)[source]#
Transforms a Gaussian Random Field into a uniform Random Field.
- Parameters:
g (array) – The values of the Gaussian Random Fields.
a (float, default=0.0) – The minimum borne of the uniform distribution.
b (float, default=1.0) – The maximum borne of the uniform distribution.
- Returns:
The uniform Random Field.
- Return type:
array
- spam.excursions.randomFields.fieldUniformToGauss(g, a=0.0, b=1.0)[source]#
Transforms a uniform Random Field into a Gaussian Random Field.
- Parameters:
g (array) – The values of the uniform Random Fields.
a (float, default=0.0) – The minimum borne of the uniform distribution.
b (float, default=1.0) – The maximum borne of the uniform distribution.
- Returns:
The Gaussian Random Field.
- Return type:
array