atmos_flux_inversion.correlations module¶
Generating correlation operators¶
-
class
atmos_flux_inversion.correlations.
HomogeneousIsotropicCorrelation
(*args, **kwargs)[source]¶ Bases:
atmos_flux_inversion.linalg.SelfAdjointLinearOperator
Homogeneous isotropic correlations using FFTs.
This embeds the physical domain passed in a larger computational domain, which allows for treatment of domains that should not be considered periodic.
Note
Do not invoke this directly. Use
HomogeneousIsotropicCorrelation.from_function()
orHomogeneousIsotropicCorrelaiton.from_array()
instead.See also
scipy.linalg.solve_circulant()
I stole the idea from here.
-
classmethod
from_array
(corr_array, is_cyclic=True)[source]¶ Create an instance with the given correlations.
- Parameters
corr_array (array_like) – The correlation of the first element of the domain with each other element.
is_cyclic (bool) –
- Returns
- Return type
-
classmethod
from_function
(corr_func, shape, is_cyclic=True)[source]¶ Create an instance to apply the correlation function.
- Parameters
corr_func (callable(dist) -> float) – The correlation of the first element of the domain with each other element.
shape (tuple of int) – The state is formally a vector, but the correlations are assumed to depend on the layout in some other shape, usually related to the physical layout. This is the other shape.
is_cyclic (bool) – Whether to assume the domain is periodic in all directions.
- Returns
- Return type
-
inv
()[source]¶ Construct the matrix inverse of this operator.
- Returns
- Return type
LinearOperator
- Raises
ValueError – if the instance is acyclic.
-
kron
(other)[source]¶ Construct the Kronecker product of this operator and other.
- Parameters
other (HomogeneousIsotropicCorrelation) – The other operator for the Kronecker product. This implementation will accept other objects, passing them along to
linalg.SchmidtKroneckerProduct
.- Returns
- Return type
-
solve
(vec)[source]¶ Solve A @ x = vec.
- Parameters
vec (array_like[N]) –
- Returns
Solution of self @ x = vec
- Return type
array_like[N]
- Raises
NotImplementedError – if the instance is acyclic.
ValueError – If vec is the wrong shape
-
atmos_flux_inversion.correlations.
make_matrix
(corr_func, shape)[source]¶ Make a correlation matrix for a domain with shape shape.
- Parameters
corr_func (callable(float) -> float[-1, 1]) – Function giving correlation between two indices a distance d from each other.
shape (tuple of int) – The underlying shape of the domain. It is viewed as a vector here, but may be more naturally seen as an N-D array. This is the shape of that array. N = prod(shape)
See also
statsmodels.stats.correlation_tools.corr_clipped()
Does something similar, and refers to other functions that may give more accurate results.
- Returns
corr – Positive definite dense array, entirely in memory
- Return type
np.ndarray[N, N]
Correlation shapes¶
Interface¶
-
class
atmos_flux_inversion.correlations.
DistanceCorrelationFunction
(length)[source]¶ Bases:
object
A correlation function that depends only on distance.
-
correlation_from_index
(*indices)[source]¶ Find the correlation between the indices.
Should be independent of the length of the underlying shape, but indices must still be even.
- Parameters
indices (tuple of int) –
- Returns
Find the correlation of the point contained in the first half of the indices with that contained in the second.
- Return type
-
Provided shapes¶
-
class
atmos_flux_inversion.correlations.
ExponentialCorrelation
(length)[source]¶ Bases:
atmos_flux_inversion.correlations.DistanceCorrelationFunction
A exponential correlation structure.
Notes
Correlation given by exp(-dist/length) where dist is the distance between the points.
-
class
atmos_flux_inversion.correlations.
BalgovindCorrelation
(length)[source]¶ Bases:
atmos_flux_inversion.correlations.DistanceCorrelationFunction
A Balgovind 3D correlation structure.
Follows Balgovind et al. 1983 recommendations for a 3D field, modified so the correlation length better matches that used by other correlation functions.
Notes
Correlation given by \((1 + 2*dist/length) exp(-2*dist/length)\)
This implementation has problems for length == 10. I have no idea why. 3 and 30 are fine.
References
Balgovind, R and Dalcher, A. and Ghil, M. and Kalnay, E. (1983). A Stochastic-Dynamic Model for the Spatial Structure of Forecast Error Statistics Monthly Weather Review 111(4) 701–722. doi:10.1175/1520-0493(1983)111<0701:asdmft>2.0.co;2
-
class
atmos_flux_inversion.correlations.
MaternCorrelation
(length, kappa=1)[source]¶ Bases:
atmos_flux_inversion.correlations.DistanceCorrelationFunction
A Matern correlation structure.
Follows Matern (1986) Spatial Variation. The specific definition seems to be similar to the parameterization in Stern’s Interpolation of Spatial Data
Notes
Correlation given by \([2^{\kappa-1}\Gamma(\kappa)]^{-1} (d/L)^{\kappa} K_{\kappa}(d/L)\) where \(\kappa\) is a smoothness parameter and \(K_{\kappa}\) is a modified Bessel function of the third kind.
References
Stein, Michael L. Interpolation of Spatial Data: Some Theory for Kridging Springer-Verlag New York. Part of Springer Series in Statistics (issn:0172-7397) isbn:978-1-4612-7166-6. doi:10.1007/978-1-4612-1494-6
-
class
atmos_flux_inversion.correlations.
GaussianCorrelation
(length)[source]¶ Bases:
atmos_flux_inversion.correlations.DistanceCorrelationFunction
A gaussian correlation structure.
Notes
Correlation given by exp(-dist**2 / (length**2)) where dist is the distance between the points.