Geospatial transforms
We provide a very powerful list of transforms that were designed to work seamlessly with geospatial data. These transforms are implemented in different packages depending on how they interact with geometries.
The list of supported transforms is continuously growing. The following code can be used to print an updated list in any project environment:
# packages to print type tree
using InteractiveUtils
using AbstractTrees
using TransformsBase
# packages with transforms
using GeoStats
# define the tree of types
AbstractTrees.children(T::Type) = subtypes(T)
# print all currently available transforms
AbstractTrees.print_tree(TransformsBase.Transform)
Transform
├─ GeometricTransform
│ ├─ Bridge
│ ├─ CoordinateTransform
│ │ ├─ Affine
│ │ ├─ LengthUnit
│ │ ├─ Morphological
│ │ ├─ Proj
│ │ ├─ Rotate
│ │ ├─ Scale
│ │ ├─ Stretch
│ │ └─ Translate
│ ├─ LambdaMuSmoothing
│ ├─ Repair
│ ├─ Shadow
│ ├─ Slice
│ └─ StdCoords
├─ Identity
├─ TableTransform
│ ├─ Aggregate
│ ├─ CookieCutter
│ ├─ Detrend
│ ├─ Downscale
│ ├─ FeatureTransform
│ │ ├─ EigenAnalysis
│ │ ├─ Remainder
│ │ ├─ ColwiseFeatureTransform
│ │ │ ├─ Center
│ │ │ ├─ Coalesce
│ │ │ ├─ LowHigh
│ │ │ ├─ Quantile
│ │ │ └─ ZScore
│ │ └─ StatelessFeatureTransform
│ │ ├─ AbsoluteUnits
│ │ ├─ Assert
│ │ ├─ Closure
│ │ ├─ Coerce
│ │ ├─ ColTable
│ │ ├─ Compose
│ │ ├─ DropConstant
│ │ ├─ DropExtrema
│ │ ├─ DropMissing
│ │ ├─ DropNaN
│ │ ├─ DropUnits
│ │ ├─ Filter
│ │ ├─ Functional
│ │ ├─ Indicator
│ │ ├─ Learn
│ │ ├─ Levels
│ │ ├─ Map
│ │ ├─ OneHot
│ │ ├─ ProjectionPursuit
│ │ ├─ Reject
│ │ ├─ Rename
│ │ ├─ Replace
│ │ ├─ RowTable
│ │ ├─ Sample
│ │ ├─ Satisfies
│ │ ├─ Select
│ │ ├─ Sort
│ │ ├─ StdFeats
│ │ ├─ StdNames
│ │ ├─ LogRatio
│ │ │ ├─ ALR
│ │ │ ├─ CLR
│ │ │ └─ ILR
│ │ ├─ Unit
│ │ └─ Unitify
│ ├─ ClusteringTransform
│ │ ├─ GHC
│ │ ├─ GSC
│ │ └─ SLIC
│ ├─ Interpolate
│ ├─ InterpolateMissing
│ ├─ InterpolateNaN
│ ├─ InterpolateNeighbors
│ ├─ Potrace
│ ├─ Rasterize
│ ├─ Simulate
│ ├─ ParallelTableTransform
│ ├─ Transfer
│ ├─ UniqueCoords
│ └─ Upscale
└─ SequentialTransform
Transforms at the leaves of the tree above should have a docstring with more information on the available options. For example, the documentation of the Select
transform is shown below:
TableTransforms.Select
— TypeSelect(col₁, col₂, ..., colₙ)
Select([col₁, col₂, ..., colₙ])
Select((col₁, col₂, ..., colₙ))
The transform that selects columns col₁
, col₂
, ..., colₙ
.
Select(col₁ => newcol₁, col₂ => newcol₂, ..., colₙ => newcolₙ)
Selects the columns col₁
, col₂
, ..., colₙ
and rename them to newcol₁
, newcol₂
, ..., newcolₙ
.
Select(regex)
Selects the columns that match with regex
.
Examples
Select(1, 3, 5)
Select([:a, :c, :e])
Select(("a", "c", "e"))
Select(1 => :x, 3 => :y)
Select(:a => :x, :b => :y)
Select("a" => "x", "b" => "y")
Select(r"[ace]")
Transforms of type FeatureTransform
operate on the attribute table, whereas transforms of type GeometricTransform
operate on the underlying geospatial domain:
TableTransforms.FeatureTransform
— TypeFeatureTransform
A transform that operates on the columns of the table containing features, i.e., simple attributes such as numbers, strings, etc.
Meshes.GeometricTransform
— TypeGeometricTransform
A method to transform the geometry (e.g. coordinates) of objects. See https://en.wikipedia.org/wiki/Geometric_transformation.
Other transforms such as Detrend
are defined in terms of both the geospatial domain and the attribute table. All transforms and pipelines implement the following functions:
TransformsBase.isrevertible
— Functionisrevertible(transform)
Tells whether or not the transform
is revertible, i.e. supports a revert
function. Defaults to false
for new transform types.
Transforms can be revertible and yet don't be invertible. Invertibility is a mathematical concept, whereas revertibility is a computational concept.
See also isinvertible
.
TransformsBase.isinvertible
— Functionisinvertible(transform)
Tells whether or not the transform
is invertible, i.e. whether it implements the inverse
function. Defaults to false
for new transform types.
Transforms can be invertible in the mathematical sense, i.e., there exists a one-to-one mapping between input and output spaces.
See also inverse
, isrevertible
.
TransformsBase.inverse
— FunctionTransformsBase.apply
— Functionnewobject, cache = apply(transform, object)
Apply transform
on the object
. Return the newobject
and a cache
to revert the transform later.
TransformsBase.revert
— Functionobject = revert(transform, newobject, cache)
Revert the transform
on the newobject
using the cache
from the corresponding apply
call and return the original object
. Only defined when the transform
isrevertible
.
TransformsBase.reapply
— Functionnewobject = reapply(transform, object, cache)
Reapply the transform
to (a possibly different) object
using a cache
that was created with a previous apply
call. Fallback to apply
without using the cache
.
Feature transforms
Please check the TableTransforms.jl documentation for an updated list of feature transforms. As an example consider the following features over a Cartesian grid and their statistics:
using DataFrames
# table of features and domain
tab = DataFrame(a=rand(1000), b=randn(1000), c=rand(1000))
dom = CartesianGrid(100, 100)
# georeference table onto domain
Ω = georef(tab, dom)
# describe features
describe(values(Ω))
Row | variable | mean | min | median | max | nmissing | eltype |
---|---|---|---|---|---|---|---|
Symbol | Float64 | Float64 | Float64 | Float64 | Int64 | DataType | |
1 | a | 0.485276 | 0.00120786 | 0.484707 | 0.999245 | 0 | Float64 |
2 | b | 0.0176609 | -2.95863 | 0.0136271 | 3.62072 | 0 | Float64 |
3 | c | 0.50863 | 7.30319e-5 | 0.511608 | 0.999331 | 0 | Float64 |
We can create a pipeline that transforms the features to their normal quantile (or scores):
pipe = Quantile()
Ω̄, cache = apply(pipe, Ω)
describe(values(Ω̄))
Row | variable | mean | min | median | max | nmissing | eltype |
---|---|---|---|---|---|---|---|
Symbol | Float64 | Float64 | Float64 | Float64 | Int64 | DataType | |
1 | a | 0.00309023 | -3.09023 | 0.00125332 | 3.09023 | 0 | Float64 |
2 | b | 0.00309023 | -3.09023 | 0.00125332 | 3.09023 | 0 | Float64 |
3 | c | 0.00309023 | -3.09023 | 0.00125332 | 3.09023 | 0 | Float64 |
We can then revert the transform given any new geospatial data in the transformed sample space:
Ωₒ = revert(pipe, Ω̄, cache)
describe(values(Ωₒ))
Row | variable | mean | min | median | max | nmissing | eltype |
---|---|---|---|---|---|---|---|
Symbol | Float64 | Float64 | Float64 | Float64 | Int64 | DataType | |
1 | a | 0.485759 | 0.00226717 | 0.48481 | 0.998913 | 0 | Float64 |
2 | b | 0.0206166 | -2.77883 | 0.0140045 | 3.60016 | 0 | Float64 |
3 | c | 0.509138 | 0.00012305 | 0.512828 | 0.998951 | 0 | Float64 |
The Learn
transform is another important transform from StatsLearnModels.jl:
StatsLearnModels.Learn
— TypeLearn(train, model, incols => outcols)
Fits the statistical learning model
using the input columns, selected by incols
, and the output columns, selected by outcols
, from the train
table.
The column selection can be a single column identifier (index or name), a collection of identifiers or a regular expression (regex).
Examples
Learn(train, model, [1, 2, 3] => "d")
Learn(train, model, [:a, :b, :c] => :d)
Learn(train, model, ["a", "b", "c"] => 4)
Learn(train, model, [1, 2, 3] => [:d, :e])
Learn(train, model, r"[abc]" => ["d", "e"])
For more details, consider watching our JuliaCon2021 talk:
Geometric transforms
Please check the Meshes.jl documentation for an updated list of geometric transforms. As an example consider the rotation of geospatial data over a Cartesian grid:
# geospatial domain
Ω = georef((Z=rand(10, 10),))
# apply geometric transform
Ωr = Ω |> Rotate(Angle2d(π/4))
fig = Mke.Figure(size = (800, 400))
viz(fig[1,1], Ω.geometry, color = Ω.Z)
viz(fig[1,2], Ωr.geometry, color = Ωr.Z)
fig
Geostatistical transforms
Bellow is the current list of transforms that operate on both the geometries and features of geospatial data. They are implemented in the GeoStatsBase.jl package.
UniqueCoords
GeoStatsTransforms.UniqueCoords
— TypeUniqueCoords(var₁ => agg₁, var₂ => agg₂, ..., varₙ => aggₙ)
Retain locations in data with unique coordinates.
Duplicates of a variable varᵢ
are aggregated with aggregation function aggᵢ
. If an aggregation function is not defined for variable varᵢ
, the default aggregation function will be used. Default aggregation function is mean
for continuous variables and first
otherwise.
Examples
UniqueCoords(1 => last, 2 => maximum)
UniqueCoords(:a => first, :b => minimum)
UniqueCoords("a" => last, "b" => maximum)
# point set with repeated points
p = rand(Point, 50)
Ω = georef((Z=rand(100),), [p; p])
Z | geometry |
---|---|
Continuous | Point |
[NoUnits] | 🖈 Cartesian{NoDatum} |
0.0178461 | (x: 0.992114 m, y: 0.524092 m, z: 0.301521 m) |
0.573037 | (x: 0.870208 m, y: 0.311693 m, z: 0.490476 m) |
0.429644 | (x: 0.441657 m, y: 0.392859 m, z: 0.577252 m) |
0.153999 | (x: 0.874132 m, y: 0.367076 m, z: 0.634166 m) |
0.242723 | (x: 0.459546 m, y: 0.0214615 m, z: 0.508143 m) |
0.423039 | (x: 0.741281 m, y: 0.240406 m, z: 0.258798 m) |
0.703679 | (x: 0.602122 m, y: 0.33676 m, z: 0.282085 m) |
0.554076 | (x: 0.837868 m, y: 0.269682 m, z: 0.381218 m) |
0.977117 | (x: 0.533105 m, y: 0.943335 m, z: 0.384398 m) |
0.828324 | (x: 0.815223 m, y: 0.793588 m, z: 0.149007 m) |
⋮ | ⋮ |
# discard repeated points
𝒰 = Ω |> UniqueCoords()
Z | geometry |
---|---|
Continuous | Point |
[NoUnits] | 🖈 Cartesian{NoDatum} |
0.105255 | (x: 0.992114 m, y: 0.524092 m, z: 0.301521 m) |
0.750587 | (x: 0.870208 m, y: 0.311693 m, z: 0.490476 m) |
0.459804 | (x: 0.441657 m, y: 0.392859 m, z: 0.577252 m) |
0.180438 | (x: 0.874132 m, y: 0.367076 m, z: 0.634166 m) |
0.146086 | (x: 0.459546 m, y: 0.0214615 m, z: 0.508143 m) |
0.577357 | (x: 0.741281 m, y: 0.240406 m, z: 0.258798 m) |
0.82928 | (x: 0.602122 m, y: 0.33676 m, z: 0.282085 m) |
0.632146 | (x: 0.837868 m, y: 0.269682 m, z: 0.381218 m) |
0.622199 | (x: 0.533105 m, y: 0.943335 m, z: 0.384398 m) |
0.491146 | (x: 0.815223 m, y: 0.793588 m, z: 0.149007 m) |
⋮ | ⋮ |
Detrend
GeoStatsTransforms.Detrend
— TypeDetrend(col₁, col₂, ..., colₙ; degree=1)
Detrend([col₁, col₂, ..., colₙ]; degree=1)
Detrend((col₁, col₂, ..., colₙ); degree=1)
The transform that detrends columns col₁
, col₂
, ..., colₙ
with a polynomial of given degree
.
Detrend(regex; degree=1)
Detrends the columns that match with regex
.
Examples
Detrend(1, 3, 5)
Detrend([:a, :c, :e])
Detrend(("a", "c", "e"))
Detrend(r"[ace]", degree=2)
Detrend(:)
References
- Menafoglio, A., Secchi, P. 2013. A Universal Kriging predictor for spatially dependent functional data of a Hilbert Space
# quadratic trend + random noise
r = range(-1, stop=1, length=100)
μ = [x^2 + y^2 for x in r, y in r]
ϵ = 0.1rand(100, 100)
Ω = georef((Z=μ+ϵ,))
# detrend and obtain noise component
𝒩 = Ω |> Detrend(:Z, degree=2)
fig = Mke.Figure(size = (800, 400))
viz(fig[1,1], Ω.geometry, color = Ω.Z)
viz(fig[1,2], 𝒩.geometry, color = 𝒩.Z)
fig
Potrace
GeoStatsTransforms.Potrace
— TypePotrace(mask; [ϵ])
Potrace(mask, var₁ => agg₁, ..., varₙ => aggₙ; [ϵ])
Trace polygons on 2D image data with Selinger's Potrace algorithm.
The categories stored in column mask
are converted into binary masks, which are then traced into multi-polygons. When provided, the option ϵ
is forwarded to Selinger's simplification algorithm.
Duplicates of a variable varᵢ
are aggregated with aggregation function aggᵢ
. If an aggregation function is not defined for variable varᵢ
, the default aggregation function will be used. Default aggregation function is mean
for continuous variables and first
otherwise.
Examples
Potrace(:mask, ϵ=0.1)
Potrace(1, 1 => last, 2 => maximum)
Potrace(:mask, :a => first, :b => minimum)
Potrace("mask", "a" => last, "b" => maximum)
References
- Selinger, P. 2003. Potrace: A polygon-based tracing algorithm
# continuous feature
Z = [sin(i/10) + sin(j/10) for i in 1:100, j in 1:100]
# binary mask
M = Z .> 0
# georeference data
Ω = georef((Z=Z, M=M))
# trace polygons using mask
𝒯 = Ω |> Potrace(:M)
fig = Mke.Figure(size = (800, 400))
viz(fig[1,1], Ω.geometry, color = Ω.Z)
viz(fig[1,2], 𝒯.geometry, color = 𝒯.Z)
fig
𝒯.geometry
2 GeometrySet
├─ Multi(4×PolyArea)
└─ Multi(4×PolyArea)
Rasterize
GeoStatsTransforms.Rasterize
— TypeRasterize(grid)
Rasterize(grid, var₁ => agg₁, ..., varₙ => aggₙ)
Rasterize geometries within specified grid
.
Rasterize(nx, ny)
Rasterize(nx, ny, var₁ => agg₁, ..., varₙ => aggₙ)
Alternatively, use the grid with size nx
by ny
obtained with discretization of the bounding box.
Duplicates of a variable varᵢ
are aggregated with aggregation function aggᵢ
. If an aggregation function is not defined for variable varᵢ
, the default aggregation function will be used. Default aggregation function is mean
for continuous variables and first
otherwise.
Examples
grid = CartesianGrid(10, 10)
Rasterize(grid)
Rasterize(10, 10)
Rasterize(grid, 1 => last, 2 => maximum)
Rasterize(10, 10, 1 => last, 2 => maximum)
Rasterize(grid, :a => first, :b => minimum)
Rasterize(10, 10, :a => first, :b => minimum)
Rasterize(grid, "a" => last, "b" => maximum)
Rasterize(10, 10, "a" => last, "b" => maximum)
A = [1, 2, 3, 4, 5]
B = [1.1, 2.2, 3.3, 4.4, 5.5]
p1 = PolyArea((2, 0), (6, 2), (2, 2))
p2 = PolyArea((0, 6), (3, 8), (0, 10))
p3 = PolyArea((3, 6), (9, 6), (9, 9), (6, 9))
p4 = PolyArea((7, 0), (10, 0), (10, 4), (7, 4))
p5 = PolyArea((1, 3), (5, 3), (6, 6), (3, 8), (0, 6))
gt = georef((; A, B), [p1, p2, p3, p4, p5])
nt = gt |> Rasterize(20, 20)
viz(nt.geometry, color = nt.A)
Clustering
Unlike traditional clustering algorithms in machine learning, geostatistical clustering (a.k.a. domaining) algorithms consider both the features and the geospatial coordinates of the data.
Consider the following data as an example:
Ω = georef((Z=[10sin(i/10) + j for i in 1:4:100, j in 1:4:100],))
viz(Ω.geometry, color = Ω.Z)
GeoStatsTransforms.GHC
— TypeGHC(k, λ; kern=:epanechnikov, link=:ward, as=:cluster)
A transform for partitioning geospatial data into k
clusters according to a range λ
using Geostatistical Hierarchical Clustering (GHC). The larger the range the more connected are nearby samples.
Parameters
k
- Approximate number of clustersλ
- Approximate range of kernel function in length unitskern
- Kernel function (:uniform
,:triangular
or:epanechnikov
)link
- Linkage function (:single
,:average
,:complete
,:ward
or:ward_presquared
)as
- Variable name used to store clustering results
References
- Fouedjio, F. 2016. A hierarchical clustering method for multivariate geostatistical data
Notes
- The range parameter controls the sparsity pattern of the pairwise distances, which can greatly affect the computational performance of the GHC algorithm. We recommend choosing a range that is small enough to connect nearby samples. For example, clustering data over a 100x100 Cartesian grid with unit spacing is possible with
λ=1.0
orλ=2.0
but the problem starts to become computationally unfeasible aroundλ=10.0
due to the density of points.
𝒞 = Ω |> GHC(20, 1.0)
viz(𝒞.geometry, color = 𝒞.cluster)
GeoStatsTransforms.GSC
— TypeGSC(k, m; σ=1.0, tol=1e-4, maxiter=10, weights=nothing, as=:cluster)
A transform for partitioning geospatial data into k
clusters using Geostatistical Spectral Clustering (GSC).
Parameters
k
- Desired number of clustersm
- Multiplicative factor for adjacent weightsσ
- Standard deviation for exponential model (default to1.0
)tol
- Tolerance of k-means algorithm (default to1e-4
)maxiter
- Maximum number of iterations (default to10
)weights
- Dictionary with weights for each attribute (default tonothing
)as
- Variable name used to store clustering results
References
- Romary et al. 2015. Unsupervised classification of multivariate geostatistical data: Two algorithms
Notes
- The algorithm implemented here is slightly different than the algorithm
described in Romary et al. 2015. Instead of setting Wᵢⱼ = 0 when i <-/-> j, we simply magnify the weight by a multiplicative factor Wᵢⱼ *= m when i <–> j. This leads to dense matrices but also better results in practice.
𝒞 = Ω |> GSC(50, 2.0)
viz(𝒞.geometry, color = 𝒞.cluster)
GeoStatsTransforms.SLIC
— TypeSLIC(k, m; tol=1e-4, maxiter=10, weights=nothing, as=:cluster)
A transform for clustering geospatial data into approximately k
clusters using Simple Linear Iterative Clustering (SLIC). The transform produces clusters of samples that are spatially connected based on a distance dₛ
and that, at the same time, are similar in terms of vars
with distance dᵥ
. The tradeoff is controlled with a hyperparameter parameter m
in an additive model dₜ = √(dᵥ² + m²(dₛ/s)²)
.
Parameters
k
- Approximate number of clustersm
- Hyperparameter of SLIC modeltol
- Tolerance of k-means algorithm (default to1e-4
)maxiter
- Maximum number of iterations (default to10
)weights
- Dictionary with weights for each attribute (default tonothing
)as
- Variable name used to store clustering results
References
- Achanta et al. 2011. SLIC superpixels compared to state-of-the-art superpixel methods
𝒞 = Ω |> SLIC(50, 0.01)
viz(𝒞.geometry, color = 𝒞.cluster)
Interpolate
GeoStatsTransforms.Interpolate
— TypeInterpolate(domain, vars₁ => model₁, ..., varsₙ => modelₙ; [parameters])
Interpolate([g₁, g₂, ..., gₙ], vars₁ => model₁, ..., varsₙ => modelₙ; [parameters])
Interpolate geospatial data on given domain
or vector of geometries [g₁, g₂, ..., gₙ]
, using geostatistical models model₁
, ..., modelₙ
for variables vars₁
, ..., varsₙ
.
Interpolate(domain, model=NN(); [parameters])
Interpolate([g₁, g₂, ..., gₙ], model=NN(); [parameters])
Interpolate geospatial data on given domain
or vector of geometries [g₁, g₂, ..., gₙ]
, using geostatistical model
for all variables.
Parameters
point
- Perform interpolation on point support (default totrue
)prob
- Perform probabilistic interpolation (default tofalse
)
See also InterpolateNeighbors
, InterpolateMissing
, InterpolateNaN
.
GeoStatsTransforms.InterpolateNeighbors
— TypeInterpolateNeighbors(domain, vars₁ => model₁, ..., varsₙ => modelₙ; [parameters])
InterpolateNeighbors([g₁, g₂, ..., gₙ], vars₁ => model₁, ..., varsₙ => modelₙ; [parameters])
Interpolate geospatial data on given domain
or set of geometries g₁
, g₂
, ..., gₙ
, using geostatistical models model₁
, ..., modelₙ
for variables vars₁
, ..., varsₙ
.
InterpolateNeighbors(domain, model=NN(); [parameters])
InterpolateNeighbors([g₁, g₂, ..., gₙ], model=NN(); [parameters])
Interpolate geospatial data on given domain
or set of geometries g₁
, g₂
, ..., gₙ
, using geostatistical model
for all variables.
Unlike Interpolate
, this transform uses neighbor search methods to fit geostatistical models at each interpolation location with a reduced number of measurements.
Parameters
minneighbors
- Minimum number of neighbors (default to1
)maxneighbors
- Maximum number of neighbors (default to10
)neighborhood
- Search neighborhood (default tonothing
)distance
- A distance defined in Distances.jl (default toEuclidean()
)point
- Perform interpolation on point support (default totrue
)prob
- Perform probabilistic interpolation (default tofalse
)
The maxneighbors
parameter can be used to perform interpolation with a subset of measurements per prediction location. If maxneighbors
is not provided, then all measurements are used.
Two neighborhood
search methods are available:
If a
neighborhood
is provided, local prediction is performed by sliding theneighborhood
in the domain.If a
neighborhood
is not provided, the prediction is performed usingmaxneighbors
nearest neighbors according todistance
.
See also Interpolate
, InterpolateMissing
, InterpolateNaN
.
GeoStatsTransforms.InterpolateMissing
— TypeInterpolateMissing(vars₁ => model₁, ..., varsₙ => modelₙ; [parameters])
InterpolateMissing(vars₁ => model₁, ..., varsₙ => modelₙ; [parameters])
Interpolate geospatial data on its own domain, using geostatistical models model₁
, ..., modelₙ
and non-missing values of the variables vars₁
, ..., varsₙ
.
InterpolateMissing(model=NN(); [parameters])
InterpolateMissing(model=NN(); [parameters])
Interpolate geospatial data on its own domain, using geostatistical model
and non-missing values of all variables.
Just like InterpolateNeighbors
, this transform uses neighbor search methods to fit geostatistical models at each interpolation location with a reduced number of measurements.
Parameters
minneighbors
- Minimum number of neighbors (default to1
)maxneighbors
- Maximum number of neighbors (default to10
)neighborhood
- Search neighborhood (default tonothing
)distance
- A distance defined in Distances.jl (default toEuclidean()
)point
- Perform interpolation on point support (default totrue
)prob
- Perform probabilistic interpolation (default tofalse
)
The maxneighbors
parameter can be used to perform interpolation with a subset of measurements per prediction location. If maxneighbors
is not provided, then all measurements are used.
Two neighborhood
search methods are available:
If a
neighborhood
is provided, local prediction is performed by sliding theneighborhood
in the domain.If a
neighborhood
is not provided, the prediction is performed usingmaxneighbors
nearest neighbors according todistance
.
See also InterpolateNaN
, InterpolateNeighbors
, Interpolate
.
GeoStatsTransforms.InterpolateNaN
— TypeInterpolateNaN(vars₁ => model₁, ..., varsₙ => modelₙ; [parameters])
InterpolateNaN(vars₁ => model₁, ..., varsₙ => modelₙ; [parameters])
Interpolate geospatial data on its own domain, using geostatistical models model₁
, ..., modelₙ
and non-NaN values of the variables vars₁
, ..., varsₙ
.
InterpolateNaN(model=NN(); [parameters])
InterpolateNaN(model=NN(); [parameters])
Interpolate geospatial data on its own domain, using geostatistical model
and non-NaN values of all variables.
Just like InterpolateNeighbors
, this transform uses neighbor search methods to fit geostatistical models at each interpolation location with a reduced number of measurements.
Parameters
minneighbors
- Minimum number of neighbors (default to1
)maxneighbors
- Maximum number of neighbors (default to10
)neighborhood
- Search neighborhood (default tonothing
)distance
- A distance defined in Distances.jl (default toEuclidean()
)point
- Perform interpolation on point support (default totrue
)prob
- Perform probabilistic interpolation (default tofalse
)
The maxneighbors
parameter can be used to perform interpolation with a subset of measurements per prediction location. If maxneighbors
is not provided, then all measurements are used.
Two neighborhood
search methods are available:
If a
neighborhood
is provided, local prediction is performed by sliding theneighborhood
in the domain.If a
neighborhood
is not provided, the prediction is performed usingmaxneighbors
nearest neighbors according todistance
.
See also InterpolateMissing
, InterpolateNeighbors
, Interpolate
.
table = (; Z=[1.,0.,1.])
coord = [(25.,25.), (50.,75.), (75.,50.)]
geotable = georef(table, coord)
grid = CartesianGrid(100, 100)
model = Kriging(GaussianVariogram(range=35.))
interp = geotable |> Interpolate(grid, model)
viz(interp.geometry, color = interp.Z)
Simulate
GeoStatsTransforms.Simulate
— TypeSimulate(domain, vars₁ => process₁, ..., varsₙ => processₙ; [parameters])
Simulate(domain, nreals, vars₁ => process₁, ..., varsₙ => processₙ; [parameters])
Simulate([g₁, g₂, ..., gₙ], vars₁ => process₁, ..., varsₙ => processₙ; [parameters])
Simulate([g₁, g₂, ..., gₙ], nreals, vars₁ => process₁, ..., varsₙ => processₙ; [parameters])
Simulate nreals
realizations of variables varsᵢ
with geostatistical process processᵢ
over given domain
or vector of geometries [g₁, g₂, ..., gₙ]
.
The parameters
are forwarded to the rand
method of the geostatistical processes.
CookieCutter
GeoStatsTransforms.CookieCutter
— TypeCookieCutter(domain, parent => process, var₁ => procmap₁, ..., varₙ => procmapₙ; [parameters])
CookieCutter(domain, nreals, parent => process, var₁ => procmap₁, ..., varₙ => procmapₙ; [parameters])
Simulate nreals
realizations of variable parent
with geostatistical process process
, and each child variable varsᵢ
with process map procmapᵢ
, over given domain
.
The process map must be an iterable of pairs of the form: value => process. Each process in the map is related to a value of the parent
realization, therefore the values of the child variables will be chosen according to the values of the corresponding parent
realization.
The parameters
are forwarded to the rand
method of the geostatistical processes.
Examples
parent = QuiltingProcess(trainimg, (30, 30))
child0 = GaussianProcess(SphericalVariogram(range=20.0, sill=0.2))
child1 = GaussianProcess(SphericalVariogram(MetricBall((200.0, 20.0))))
transform = CookieCutter(domain, :parent => parent, :child => [0 => child0, 1 => child1])