RingGrids
RingGrids is a submodule that has been developed for SpeedyWeather.jl which is technically independent (SpeedyWeather.jl however imports it and so does SpeedyTransforms) and can also be used without running simulations. It is just not put into its own respective repository.
RingGrids defines several iso-latitude grids, which are mathematically described in the section on Grids. In brief, they include the regular latitude-longitude grids (here called FullClenshawGrid
) as well as grids which latitudes are shifted to the Gaussian latitudes and reduced grids, meaning that they have a decreasing number of longitudinal points towards the poles to be more equal-area than full grids.
RingGrids defines and exports the following grids:
- full grids:
FullClenshawGrid
,FullGaussianGrid
,FullHEALPix
, andFullOctaHEALPix
- reduced grids:
OctahedralGaussianGrid
,OctahedralClenshawGrid
,OctaHEALPixGrid
andHEALPixGrid
The following explanation of how to use these can be mostly applied to any of them, however, there are certain functions that are not defined, e.g. the full grids can be trivially converted to a Matrix
(i.e. they are rectangular grids) but not the OctahedralGaussianGrid
.
We use the term ring, short for iso-latitude ring, to refer to a sequence of grid points that all share the same latitude. A latitude-longitude grid is a ring grid, as it organises its grid-points into rings. However, other grids, like the cubed-sphere are not based on iso-latitude rings. SpeedyWeather.jl only works with ring grids because its a requirement for the Spherical Harmonic Transform to be efficient. See Grids.
Creating data on a RingGrid
Every grid
in RingGrids has a grid.data
field, which is a vector containing the data on the grid. The grid points are unravelled west to east then north to south, meaning that it starts at 90˚N and 0˚E then walks eastward for 360˚ before jumping on the next latitude ring further south, this way circling around the sphere till reaching the south pole. This may also be called ring order.
Data in a Matrix
which follows this ring order can be put on a FullGaussianGrid
like so
using SpeedyWeather.RingGrids
map = randn(Float32, 8, 4)
8×4 Matrix{Float32}:
-0.623103 -1.60846 -0.261529 -0.994688
0.914091 -0.428023 -0.648645 -3.07402
-1.68958 0.77597 0.842834 0.0387439
-0.648679 -0.720717 -0.539725 -1.81425
0.166466 0.993667 0.753921 0.708371
-0.705848 -0.054411 2.25518 0.658927
-2.29844 -1.01451 -1.24396 -0.839357
-0.124222 -0.505617 -1.43198 2.36439
grid = FullGaussianGrid(map, input_as=Matrix)
32-element, 4-ring FullGaussianGrid{Float32}:
-0.6231032
0.9140912
-1.6895838
-0.64867884
0.16646568
-0.7058484
-2.2984436
-0.12422195
-1.6084567
-0.42802262
⋮
-1.4319807
-0.9946878
-3.0740242
0.038743917
-1.8142525
0.7083713
0.65892684
-0.8393567
2.3643882
Note that input_as=Matrix
is necessary as, RingGrids have a flattened horizontal dimension into a vector. To distinguish the 2nd horizontal dimension from a possible vertical dimension the keyword argument here is required.
A full Gaussian grid has always $2N$ x $N$ grid points, but a FullClenshawGrid
has $2N$ x $N-1$, if those dimensions don't match, the creation will throw an error. To reobtain the data from a grid, you can access its data
field which returns a normal Vector
grid.data
32-element Vector{Float32}:
-0.6231032
0.9140912
-1.6895838
-0.64867884
0.16646568
-0.7058484
-2.2984436
-0.12422195
-1.6084567
-0.42802262
⋮
-1.4319807
-0.9946878
-3.0740242
0.038743917
-1.8142525
0.7083713
0.65892684
-0.8393567
2.3643882
Which can be reshaped to reobtain map
from above. Alternatively you can Matrix(grid)
to do this in one step
map == Matrix(FullGaussianGrid(map))
false
You can also use zeros
, ones
, rand
, randn
to create a grid, whereby nlat_half
, i.e. the number of latitude rings on one hemisphere, Equator included, is used as a resolution parameter and here as a second argument.
nlat_half = 4
grid = randn(OctahedralGaussianGrid{Float16}, nlat_half)
208-element, 8-ring OctahedralGaussianGrid{Float16}:
0.7275
-0.129
0.831
-0.4302
0.719
0.4917
-1.281
-1.081
1.618
0.4082
⋮
-1.056
0.2205
-0.857
-0.7534
-0.1449
-0.2091
-1.894
0.04932
-0.9897
and any element type T
can be used for OctahedralGaussianGrid{T}
and similar for other grid types.
Visualising RingGrid data
As only the full grids can be reshaped into a matrix, the underlying data structure of any AbstractGrid
is a vector. As shown in the examples above, one can therefore inspect the data as if it was a vector. But as that data has, through its <:AbstractGrid
type, all the geometric information available to plot it on a map, RingGrids also exports plot
function, based on UnicodePlots' heatmap
.
nlat_half = 24
grid = randn(OctahedralGaussianGrid, nlat_half)
RingGrids.plot(grid)
48-ring OctahedralGaussianGrid{Float64}
┌────────────────────────────────────────────────────────────┐ 3
90 │▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ ┌──┐
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
˚N │▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
│▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ │▄▄│
-90 │▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄▄│ └──┘
└────────────────────────────────────────────────────────────┘ -3
0 ˚E 360
(Note that to skip the RingGrids.
in the last line you can do import SpeedyWeather.RingGrids: plot
, import SpeedyWeather: plot
or simply using SpeedyWeather
. It's just that LowerTriangularMatrices
also defines plot
which otherwise causes naming conflicts.)
Indexing RingGrids
All RingGrids have a single index ij
which follows the ring order. While this is obviously not super exciting here are some examples how to make better use of the information that the data sits on a grid.
We obtain the latitudes of the rings of a grid by calling get_latd
(get_lond
is only defined for full grids, or use get_latdlonds
for latitudes, longitudes per grid point not per ring)
grid = randn(OctahedralClenshawGrid, 5)
latd = get_latd(grid)
9-element Vector{Float64}:
72.0
54.0
36.0
18.0
0.0
-18.0
-36.0
-54.0
-72.0
Now we could calculate Coriolis and add it on the grid as follows
rotation = 7.29e-5 # angular frequency of Earth's rotation [rad/s]
coriolis = 2rotation*sind.(latd) # vector of coriolis parameters per latitude ring
rings = eachring(grid)
for (j, ring) in enumerate(rings)
f = coriolis[j]
for ij in ring
grid[ij] += f
end
end
eachring
creates a vector of UnitRange
indices, such that we can loop over the ring index j
(j=1
being closest to the North pole) pull the coriolis parameter at that latitude and then loop over all in-ring indices i
(changing longitudes) to do something on the grid. Something similar can be done to scale/unscale with the cosine of latitude for example. We can always loop over all grid-points like so
for ij in eachgridpoint(grid)
grid[ij]
end
or use eachindex
instead.
Interpolation on RingGrids
In most cases we will want to use RingGrids so that our data directly comes with the geometric information of where the grid-point is one the sphere. We have seen how to use get_latd
, get_lond
, ... for that above. This information generally can also be used to interpolate our data from grid to another or to request an interpolated value on some coordinates. Using our data on grid
which is an OctahedralGaussianGrid
from above we can use the interpolate
function to get it onto a FullGaussianGrid
(or any other grid for purpose)
grid = randn(OctahedralGaussianGrid{Float32}, 4)
208-element, 8-ring OctahedralGaussianGrid{Float32}:
0.376631
-0.89818704
-0.34631652
-0.15891403
0.5113424
-0.1246499
-1.6417928
0.47886133
-0.63897145
1.1653304
⋮
0.64781415
-0.062504895
-1.1181656
-0.4209753
-0.12772472
-0.47776118
0.56537396
0.8373177
-1.1051904
interpolate(FullGaussianGrid, grid)
128-element, 8-ring FullGaussianGrid{Float64}:
0.3766309916973114
-0.7602194100618362
-0.2526152729988098
0.3437782973051071
-0.12464989721775055
-1.1116292476654053
-0.08005505800247192
0.7142549455165863
1.9196585416793823
0.2073495090007782
⋮
1.376798689365387
1.09537935256958
0.47023439407348633
-0.5903352573513985
-0.5952728986740125
-0.12772472202777863
-0.21697738766670227
0.701345831155777
-0.6195633709430695
By default this will linearly interpolate (it's an Anvil interpolator, see below) onto a grid with the same nlat_half
, but we can also coarse-grain or fine-grain by specifying nlat_half
directly as 2nd argument
interpolate(FullGaussianGrid, 6, grid)
288-element, 12-ring FullGaussianGrid{Float64}:
0.290863260175366
-0.43175860280211065
-0.3260237418584422
-0.13715744772762584
0.07855134261370679
0.310393976189551
-0.05011389450870915
-0.910095922488001
-0.12042974622985006
-0.01977996052983763
⋮
-0.3533210665116874
-0.5542784703679915
-0.2048750173775651
-0.038647941025245214
-0.2370636941549111
0.19628834853405822
0.5252955974342512
0.1773463146938694
-0.39448283697645725
So we got from an 8-ring OctahedralGaussianGrid{Float16}
to a 12-ring FullGaussianGrid{Float64}
, so it did a conversion from Float16
to Float64
on the fly too, because the default precision is Float64
unless specified. interpolate(FullGaussianGrid{Float16}, 6, grid)
would have interpolated onto a grid with element type Float16
.
One can also interpolate onto a given coordinate ˚N, ˚E like so
interpolate(30.0, 10.0, grid)
0.10765388f0
we interpolated the data from grid
onto 30˚N, 10˚E. To do this simultaneously for many coordinates they can be packed into a vector too
interpolate([30.0, 40.0, 50.0], [10.0, 10.0, 10.0], grid)
3-element Vector{Float32}:
0.10765388
0.25872988
0.4027013
which returns the data on grid
at 30˚N, 40˚N, 50˚N, and 10˚E respectively. Note how the interpolation here retains the element type of grid
.
Performance for RingGrid interpolation
Every time an interpolation like interpolate(30.0, 10.0, grid)
is called, several things happen, which are important to understand to know how to get the fastest interpolation out of this module in a given situation. Under the hood an interpolation takes three arguments
- output vector
- input grid
- interpolator
The output vector is just an array into which the interpolated data is written, providing this prevents unnecessary allocation of memory in case the destination array of the interpolation already exists. The input grid contains the data which is subject to interpolation, it must come on a ring grid, however, its coordinate information is actually already in the interpolator. The interpolator knows about the geometry of the grid the data is coming on and the coordinates it is supposed to interpolate onto. It has therefore precalculated the indices that are needed to access the right data on the input grid and the weights it needs to apply in the actual interpolation operation. The only thing it does not know is the actual data values of that grid. So in the case you want to interpolate from grid A to grid B many times, you can just reuse the same interpolator. If you want to change the coordinates of the output grid but its total number of points remain constants then you can update the locator inside the interpolator and only else you will need to create a new interpolator. Let's look at this in practice. Say we have two grids an want to interpolate between them
grid_in = rand(HEALPixGrid, 4)
grid_out = zeros(FullClenshawGrid, 6)
interp = RingGrids.interpolator(grid_out, grid_in)
AnvilInterpolator{Float64, HEALPixGrid}
├ from: 7-ring HEALPixGrid, 48 grid points
└ onto: 264 points
Now we have created an interpolator interp
which knows about the geometry where to interpolate from and the coordinates there to interpolate to. It is also initialized, meaning it has precomputed the indices to of grid_in
that are supposed to be used. It just does not know about the data of grid_in
(and neither of grid_out
which will be overwritten anyway). We can now do
interpolate!(grid_out, grid_in, interp)
grid_out
264-element, 11-ring FullClenshawGrid{Float64}:
0.49091025894286805
0.4187850656215237
0.34665987230017936
0.274534678978835
0.285834237068754
0.297133795158673
0.308433353248592
0.319732911338511
0.33103246942842995
0.34233202751834896
⋮
0.5721016762045251
0.5080220248633976
0.4439423735222701
0.37986272218114253
0.31578307084001506
0.25170341949888775
0.18762376815775988
0.17990336861543432
0.17218296907310876
which is identical to interpolate(grid_out, grid_in)
but you can reuse interp
for other data. The interpolation can also handle various element types (the interpolator interp
does not have to be updated for this either)
grid_out = zeros(FullClenshawGrid{Float16}, 6);
interpolate!(grid_out, grid_in, interp)
grid_out
264-element, 11-ring FullClenshawGrid{Float16}:
0.491
0.4187
0.3467
0.2744
0.286
0.297
0.3083
0.3198
0.331
0.3423
⋮
0.5723
0.508
0.4438
0.38
0.3157
0.2517
0.1876
0.1799
0.1722
and we have converted data from a HEALPixGrid{Float64}
(Float64
is always default if not specified) to a FullClenshawGrid{Float16}
including the type conversion Float64-Float16 on the fly. Technically there are three data types and their combinations possible: The input data will come with a type, the output array has an element type and the interpolator has precomputed weights with a given type. Say we want to go from Float16 data on an OctahedralGaussianGrid
to Float16 on a FullClenshawGrid
but using Float32 precision for the interpolation itself, we would do this by
grid_in = randn(OctahedralGaussianGrid{Float16}, 24)
grid_out = zeros(FullClenshawGrid{Float16}, 24)
interp = RingGrids.interpolator(Float32, grid_out, grid_in)
interpolate!(grid_out, grid_in, interp)
grid_out
4512-element, 47-ring FullClenshawGrid{Float16}:
0.2299
0.3005
0.3713
0.4421
0.5127
0.8203
0.7144
0.609
0.503
0.09174
⋮
-0.9536
-1.094
-0.984
-0.8735
-0.7637
-0.4827
-0.2325
0.0176
0.2678
As a last example we want to illustrate a situation where we would always want to interpolate onto 10 coordinates, but their locations may change. In order to avoid recreating an interpolator object we would do (AnvilInterpolator
is described in Anvil interpolator)
npoints = 10 # number of coordinates to interpolate onto
interp = AnvilInterpolator(Float32, HEALPixGrid, 24, npoints)
AnvilInterpolator{Float32, HEALPixGrid}
├ from: 47-ring HEALPixGrid, 1728 grid points
└ onto: 10 points
with the first argument being the number format used during interpolation, then the input grid type, its resolution in terms of nlat_half
and then the number of points to interpolate onto. However, interp
is not yet initialized as it does not know about the destination coordinates yet. Let's define them, but note that we already decided there's only 10 of them above.
latds = collect(0.0:5.0:45.0)
londs = collect(-10.0:2.0:8.0)
now we can update the locator inside our interpolator as follows
RingGrids.update_locator!(interp, latds, londs)
With data matching the input from above, a nlat_half=24
HEALPixGrid, and allocate 10-element output vector
output_vec = zeros(10)
grid_input = rand(HEALPixGrid, 24)
we can use the interpolator as follows
interpolate!(output_vec, grid_input, interp)
10-element Vector{Float64}:
0.9236585636846435
0.6291884083279105
0.6621649835144147
0.49904702937360523
0.7923496908950889
0.6948980871496233
0.3379956050208286
0.5517189510601717
0.540411964417737
0.7264419815083791
which is the approximately the same as doing it directly without creating an interpolator first and updating its locator
interpolate(latds, londs, grid_input)
10-element Vector{Float64}:
0.9236585647105592
0.6291884106589929
0.6621649830202012
0.4990470402930828
0.7923496893010221
0.6948980887989882
0.3379956041415177
0.5517189379166438
0.5404119774293872
0.7264419815603529
but allows for a reuse of the interpolator. Note that the two output arrays are not exactly identical because we manually set our interpolator interp
to use Float32
for the interpolation whereas the default is Float64
.
Anvil interpolator
Currently the only interpolator implemented is a 4-point bilinear interpolator, which schematically works as follows. Anvil interpolation is the bilinear average of a, b, c, d which are values at grid points in an anvil-shaped configuration at location x, which is denoted by Δab, Δcd, Δy, the fraction of distances between a-b, c-d, and ab-cd, respectively. Note that a, c and b, d do not necessarily share the same longitude/x-coordinate.
0..............1 # fraction of distance Δab between a, b
|< Δab >|
0^ a -------- o - b # anvil-shaped average of a, b, c, d at location x
.Δy |
. |
.v x
. |
1 c ------ o ---- d
|< Δcd >|
0...............1 # fraction of distance Δcd between c, d
^ fraction of distance Δy between a-b and c-d.
This interpolation is chosen as by definition of the ring grids, a and b share the same latitude, so do c and d, but the longitudes can be different for all four, a, b, c, d.
Function index
Core.Type
— Method() Initialize an instance of the grid from an Array. For keyword argument input_as=Vector
(default) the leading dimension is interpreted as a flat vector of all horizontal entries in one layer. For input_as==Matrx
the first two leading dimensions are interpreted as longitute and latitude. This is only possible for full grids that are a subtype of AbstractFullGridArray
.
SpeedyWeather.RingGrids.AbstractFullGrid
— TypeAn AbstractFullGrid
is a horizontal grid with a constant number of longitude points across latitude rings. Different latitudes can be used, Gaussian latitudes, equi-angle latitudes (also called Clenshaw from Clenshaw-Curtis quadrature), or others.
SpeedyWeather.RingGrids.AbstractFullGridArray
— TypeSubtype of AbstractGridArray
for all N-dimensional arrays of ring grids that have the same number of longitude points on every ring. As such these (horizontal) grids are representable as a matrix, with denser grid points towards the poles.
SpeedyWeather.RingGrids.AbstractGrid
— TypeAbstract supertype for all ring grids, representing 2-dimensional data on the sphere unravelled into a Julia Vector
. Subtype of AbstractGridArray
with N=1
and ArrayType=Vector{T}
of eltype T
.
SpeedyWeather.RingGrids.AbstractGridArray
— TypeAbstract supertype for all arrays of ring grids, representing N
-dimensional data on the sphere in two dimensions (but unravelled into a vector in the first dimension, the actual "ring grid") plus additional N-1
dimensions for the vertical and/or time etc. Parameter T
is the eltype
of the underlying data, held as in the array type ArrayType
(Julia's Array
for CPU or others for GPU).
Ring grids have several consecuitive grid points share the same latitude (= a ring), grid points on a given ring are equidistant. Grid points are ordered 0 to 360˚E, starting around the north pole, ring by ring to the south pole.
SpeedyWeather.RingGrids.AbstractInterpolator
— Typeabstract type AbstractInterpolator{NF, G} end
Supertype for Interpolators. Every Interpolator <: AbstractInterpolator is expected to have two fields,
- geometry, which describes the grid G to interpolate from
- locator, which locates the indices on G and their weights to interpolate onto a new grid.
NF is the number format used to calculate the interpolation, which can be different from the input data and/or the interpolated data on the new grid.
SpeedyWeather.RingGrids.AbstractLocator
— TypeAbstractLocator{NF}
Supertype of every Locator, which locates the indices on a grid to be used to perform an interpolation. E.g. AnvilLocator uses a 4-point stencil for every new coordinate to interpolate onto. Higher order stencils can be implemented by defining OtherLocator <: AbstractLocactor.
SpeedyWeather.RingGrids.AbstractReducedGrid
— TypeHorizontal abstract type for all AbstractReducedGridArray
with N=1
(i.e. horizontal only) and ArrayType
of Vector{T}
with element type T
.
SpeedyWeather.RingGrids.AbstractReducedGridArray
— TypeSubtype of AbstractGridArray
for arrays of rings grids that have a reduced number of longitude points towards the poles, i.e. they are not "full", see AbstractFullGridArray
. Data on these grids cannot be represented as matrix and has to be unravelled into a vector, ordered 0 to 360˚E then north to south, ring by ring. Examples for reduced grids are the octahedral Gaussian or Clenshaw grids, or the HEALPix grid.
SpeedyWeather.RingGrids.AbstractSphericalDistance
— Typeabstract type AbstractSphericalDistance end
Super type of formulas to calculate the spherical distance or great-circle distance. To define a NewFormula
, define struct NewFormula <: AbstractSphericalDistance end
and the actual calculation as a functor
function NewFormula(lonlat1::Tuple, lonlat2::Tuple; radius=DEFAULT_RADIUS, kwargs...)
assuming inputs in degrees and returning the distance in meters (or radians for radius=1). Then use the general interface spherical_distance(NewFormula, args...; kwargs...)
SpeedyWeather.RingGrids.AnvilLocator
— TypeAnvilLocator{NF<:AbstractFloat} <: AbtractLocator
Contains arrays that locates grid points of a given field to be uses in an interpolation and their weights. This Locator is a 4-point average in an anvil-shaped grid-point arrangement between two latitude rings.
SpeedyWeather.RingGrids.AnvilLocator
— MethodL = AnvilLocator( ::Type{NF}, # number format used for the interpolation
npoints::Integer # number of points to interpolate onto
) where {NF<:AbstractFloat}
Zero generator function for the 4-point average AnvilLocator. Use update_locator! to update the grid indices used for interpolation and their weights. The number format NF is the format used for the calculations within the interpolation, the input data and/or output data formats may differ.
SpeedyWeather.RingGrids.FullClenshawArray
— TypeA FullClenshawArray
is an array of full grid, subtyping AbstractFullGridArray
, that use equidistant latitudes for each ring (a regular lon-lat grid). These require the Clenshaw-Curtis quadrature in the spectral transform, hence the name. One ring is on the equator, total number of rings is odd, no rings on the north or south pole. First dimension of the underlying N
-dimensional data
represents the horizontal dimension, in ring order (0 to 360˚E, then north to south), other dimensions are used for the vertical and/or time or other dimensions. The resolution parameter of the horizontal grid is nlat_half
(number of latitude rings on one hemisphere, Equator included) and the ring indices are precomputed in rings
. Fields are
data::AbstractArray
nlat_half::Int64
rings::Vector{UnitRange{Int64}}
SpeedyWeather.RingGrids.FullClenshawGrid
— TypeA FullClenshawArray
but constrained to N=1
dimensions (horizontal only) and data is a Vector{T}
.
SpeedyWeather.RingGrids.FullGaussianArray
— TypeA FullGaussianArray
is an array of full grids, subtyping AbstractFullGridArray
, that use Gaussian latitudes for each ring. First dimension of the underlying N
-dimensional data
represents the horizontal dimension, in ring order (0 to 360˚E, then north to south), other dimensions are used for the vertical and/or time or other dimensions. The resolution parameter of the horizontal grid is nlat_half
(number of latitude rings on one hemisphere, Equator included) and the ring indices are precomputed in rings
. Fields are
data::AbstractArray
: Data array, west to east, ring by ring, north to south.nlat_half::Int64
: Number of latitudes on one hemisphererings::Vector{UnitRange{Int64}}
: Precomputed ring indices, ranging from first to last grid point on every ring.
SpeedyWeather.RingGrids.FullGaussianGrid
— TypeA FullGaussianArray
but constrained to N=1
dimensions (horizontal only) and data is a Vector{T}
.
SpeedyWeather.RingGrids.FullHEALPixArray
— TypeA FullHEALPixArray
is an array of full grids, subtyping AbstractFullGridArray
, that use HEALPix latitudes for each ring. This type primarily equists to interpolate data from the (reduced) HEALPixGrid onto a full grid for output.
First dimension of the underlying N
-dimensional data
represents the horizontal dimension, in ring order (0 to 360˚E, then north to south), other dimensions are used for the vertical and/or time or other dimensions. The resolution parameter of the horizontal grid is nlat_half
(number of latitude rings on one hemisphere, Equator included) and the ring indices are precomputed in rings
. Fields are
data::AbstractArray
nlat_half::Int64
rings::Vector{UnitRange{Int64}}
SpeedyWeather.RingGrids.FullHEALPixGrid
— TypeA FullHEALPixArray
but constrained to N=1
dimensions (horizontal only) and data is a Vector{T}
.
SpeedyWeather.RingGrids.FullOctaHEALPixArray
— TypeA FullOctaHEALPixArray
is an array of full grids, subtyping AbstractFullGridArray
that use OctaHEALPix latitudes for each ring. This type primarily equists to interpolate data from the (reduced) OctaHEALPixGrid onto a full grid for output.
First dimension of the underlying N
-dimensional data
represents the horizontal dimension, in ring order (0 to 360˚E, then north to south), other dimensions are used for the vertical and/or time or other dimensions. The resolution parameter of the horizontal grid is nlat_half
(number of latitude rings on one hemisphere, Equator included) and the ring indices are precomputed in rings
. Fields are
data::AbstractArray
nlat_half::Int64
rings::Vector{UnitRange{Int64}}
SpeedyWeather.RingGrids.FullOctaHEALPixGrid
— TypeA FullOctaHEALPixArray
but constrained to N=1
dimensions (horizontal only) and data is a Vector{T}
.
SpeedyWeather.RingGrids.GridGeometry
— TypeGridGeometry{G<:AbstractGrid}
contains general precomputed arrays describing the grid of G.
SpeedyWeather.RingGrids.GridGeometry
— MethodG = GridGeometry( Grid::Type{<:AbstractGrid},
nlat_half::Integer)
Precomputed arrays describing the geometry of the Grid with resolution nlat_half. Contains latitudes and longitudes of grid points, their ring index j and their unravelled indices ij.
SpeedyWeather.RingGrids.HEALPixArray
— TypeA HEALPixArray
is an array of HEALPix grids, subtyping AbstractReducedGridArray
. First dimension of the underlying N
-dimensional data
represents the horizontal dimension, in ring order (0 to 360˚E, then north to south), other dimensions are used for the vertical and/or time or other dimensions. The resolution parameter of the horizontal grid is nlat_half
(number of latitude rings on one hemisphere, Equator included) which has to be even (non-fatal error thrown otherwise) which is less strict than the original HEALPix formulation (only powers of two for nside = nlat_half/2). Ring indices are precomputed in rings
.
A HEALPix grid has 12 faces, each nside
xnside
grid points, each covering the same area of the sphere. They start with 4 longitude points on the northern-most ring, increase by 4 points per ring in the "polar cap" (the top half of the 4 northern-most faces) but have a constant number of longitude points in the equatorial belt. The southern hemisphere is symmetric to the northern, mirrored around the Equator. HEALPix grids have a ring on the Equator. For more details see Górski et al. 2005, DOI:10.1086/427976.
rings
are the precomputed ring indices, for nlat_half = 4 it is rings = [1:4, 5:12, 13:20, 21:28, 29:36, 37:44, 45:48]
. So the first ring has indices 1:4 in the unravelled first dimension, etc. For efficient looping see eachring
and eachgrid
. Fields are
data::AbstractArray
nlat_half::Int64
rings::Vector{UnitRange{Int64}}
SpeedyWeather.RingGrids.HEALPixGrid
— TypeAn HEALPixArray
but constrained to N=1
dimensions (horizontal only) and data is a Vector{T}
.
SpeedyWeather.RingGrids.Haversine
— MethodHaversine(lonlat1::Tuple, lonlat2::Tuple; radius) -> Any
Haversine formula calculating the great-circle or spherical distance (in meters) on the sphere between two tuples of longitude-latitude points in degrees ˚E, ˚N. Use keyword argument radius
to change the radius of the sphere (default 6371e3 meters, Earth's radius), use radius=1
to return the central angle in radians or radius=360/2π
to return degrees.
SpeedyWeather.RingGrids.OctaHEALPixArray
— TypeAn OctaHEALPixArray
is an array of OctaHEALPix grids, subtyping AbstractReducedGridArray
. First dimension of the underlying N
-dimensional data
represents the horizontal dimension, in ring order (0 to 360˚E, then north to south), other dimensions are used for the vertical and/or time or other dimensions. The resolution parameter of the horizontal grid is nlat_half
(number of latitude rings on one hemisphere, Equator included) and the ring indices are precomputed in rings
.
An OctaHEALPix grid has 4 faces, each nlat_half x nlat_half
in size, covering 90˚ in longitude, pole to pole. As part of the HEALPix family of grids, the grid points are equal area. They start with 4 longitude points on the northern-most ring, increase by 4 points per ring towards the Equator with one ring on the Equator before reducing the number of points again towards the south pole by 4 per ring. There is no equatorial belt for OctaHEALPix grids. The southern hemisphere is symmetric to the northern, mirrored around the Equator. OctaHEALPix grids have a ring on the Equator. For more details see Górski et al. 2005, DOI:10.1086/427976, the OctaHEALPix grid belongs to the family of HEALPix grids with Nθ = 1, Nφ = 4 but is not explicitly mentioned therein.
rings
are the precomputed ring indices, for nlat_half = 3 (in contrast to HEALPix this can be odd) it is rings = [1:4, 5:12, 13:24, 25:32, 33:36]
. For efficient looping see eachring
and eachgrid
. Fields are
data::AbstractArray
nlat_half::Int64
rings::Vector{UnitRange{Int64}}
SpeedyWeather.RingGrids.OctaHEALPixGrid
— TypeAn OctaHEALPixArray
but constrained to N=1
dimensions (horizontal only) and data is a Vector{T}
.
SpeedyWeather.RingGrids.OctahedralClenshawArray
— TypeAn OctahedralClenshawArray
is an array of octahedral grids, subtyping AbstractReducedGridArray
, that use equidistant latitudes for each ring, the same as for FullClenshawArray
. First dimension of the underlying N
-dimensional data
represents the horizontal dimension, in ring order (0 to 360˚E, then north to south), other dimensions are used for the vertical and/or time or other dimensions. The resolution parameter of the horizontal grid is nlat_half
(number of latitude rings on one hemisphere, Equator included) and the ring indices are precomputed in rings
.
These grids are called octahedral (same as for the OctahedralGaussianArray
which only uses different latitudes) because after starting with 20 points on the first ring around the north pole (default) they increase the number of longitude points for each ring by 4, such that they can be conceptually thought of as lying on the 4 faces of an octahedron on each hemisphere. Hence, these grids have 20, 24, 28, ... longitude points for ring 1, 2, 3, ... Clenshaw grids have a ring on the Equator which has 16 + 4nlat_half longitude points before reducing the number of longitude points per ring by 4 towards the southern-most ring j = nlat
. rings
are the precomputed ring indices, the the example above rings = [1:20, 21:44, 45:72, ...]
. For efficient looping see eachring
and eachgrid
. Fields are
data::AbstractArray
nlat_half::Int64
rings::Vector{UnitRange{Int64}}
SpeedyWeather.RingGrids.OctahedralClenshawGrid
— TypeAn OctahedralClenshawArray
but constrained to N=1
dimensions (horizontal only) and data is a Vector{T}
.
SpeedyWeather.RingGrids.OctahedralGaussianArray
— TypeAn OctahedralGaussianArray
is an array of octahedral grids, subtyping AbstractReducedGridArray
, that use Gaussian latitudes for each ring. First dimension of the underlying N
-dimensional data
represents the horizontal dimension, in ring order (0 to 360˚E, then north to south), other dimensions are used for the vertical and/or time or other dimensions. The resolution parameter of the horizontal grid is nlat_half
(number of latitude rings on one hemisphere, Equator included) and the ring indices are precomputed in rings
.
These grids are called octahedral because after starting with 20 points on the first ring around the north pole (default) they increase the number of longitude points for each ring by 4, such that they can be conceptually thought of as lying on the 4 faces of an octahedron on each hemisphere. Hence, these grids have 20, 24, 28, ... longitude points for ring 1, 2, 3, ... There is no ring on the Equator and the two rings around it have 16 + 4nlat_half longitude points before reducing the number of longitude points per ring by 4 towards the southern-most ring j = nlat. rings
are the precomputed ring indices, the the example above rings = [1:20, 21:44, 45:72, ...]
. For efficient looping see eachring
and eachgrid
. Fields are
data::AbstractArray
nlat_half::Int64
rings::Vector{UnitRange{Int64}}
SpeedyWeather.RingGrids.OctahedralGaussianGrid
— TypeAn OctahedralGaussianArray
but constrained to N=1
dimensions (horizontal only) and data is a Vector{T}
.
SpeedyWeather.RingGrids.OctaminimalGaussianArray
— TypeAn OctaminimalGaussianArray
is an array of octahedral grids, subtyping AbstractReducedGridArray
, that use Gaussian latitudes for each ring. First dimension of the underlying N
-dimensional data
represents the horizontal dimension, in ring order (0 to 360˚E, then north to south), other dimensions are used for the vertical and/or time or other dimensions. The resolution parameter of the horizontal grid is nlat_half
(number of latitude rings on one hemisphere, Equator included) and the ring indices are precomputed in rings
.
These grids are called octahedral because after starting with 4 points on the first ring around the north pole (default) they increase the number of longitude points for each ring by 4, such that they can be conceptually thought of as lying on the 4 faces of an octahedron on each hemisphere. Hence, these grids have 4, 8, 12, ... longitude points for ring 1, 2, 3, ... which is in contrast to the OctahedralGaussianArray
which starts with 20 points around the poles, hence "minimal". There is no ring on the Equator and the two rings around it have 4nlat_half longitude points before reducing the number of longitude points per ring by 4 towards the southern-most ring j = nlat. rings
are the precomputed ring indices, in the example above rings = [1:4, 5:12, 13:24, ...]
. For efficient looping see eachring
and eachgrid
. Fields are
data::AbstractArray
nlat_half::Int64
rings::Vector{UnitRange{Int64}}
SpeedyWeather.RingGrids.OctaminimalGaussianGrid
— TypeAn OctaminimalGaussianArray
but constrained to N=1
dimensions (horizontal only) and data is a Vector{T}
.
Base.sizeof
— Methodsizeof(G::AbstractGridArray) -> Any
Size of underlying data array plus precomputed ring indices.
SpeedyWeather.RingGrids.Matrix!
— MethodMatrix!(M::AbstractMatrix,
G::OctaHEALPixGrid;
quadrant_rotation=(0, 1, 2, 3),
matrix_quadrant=((2, 2), (1, 2), (1, 1), (2, 1)),
)
Sorts the gridpoints in G
into the matrix M
without interpolation. Every quadrant of the grid G
is rotated as specified in quadrant_rotation
, 0 is no rotation, 1 is 90˚ clockwise, 2 is 180˚ etc. Grid quadrants are counted eastward starting from 0˚E. The grid quadrants are moved into the matrix quadrant (i, j) as specified. Defaults are equivalent to centered at 0˚E and a rotation such that the North Pole is at M's midpoint.
SpeedyWeather.RingGrids.Matrix!
— MethodMatrix!(
M::AbstractMatrix,
G::OctahedralClenshawGrid;
kwargs...
) -> AbstractMatrix
Sorts the gridpoints in G
into the matrix M
without interpolation. Every quadrant of the grid G
is rotated as specified in quadrant_rotation
, 0 is no rotation, 1 is 90˚ clockwise, 2 is 180˚ etc. Grid quadrants are counted eastward starting from 0˚E. The grid quadrants are moved into the matrix quadrant (i, j) as specified. Defaults are equivalent to centered at 0˚E and a rotation such that the North Pole is at M's midpoint.
SpeedyWeather.RingGrids.Matrix!
— MethodMatrix!(MGs::Tuple{AbstractMatrix{T}, OctaHEALPixGrid}...; kwargs...)
Like Matrix!(::AbstractMatrix, ::OctaHEALPixGrid)
but for simultaneous processing of tuples ((M1, G1), (M2, G2), ...)
with matrices Mi
and grids Gi
. All matrices and grids have to be of the same size respectively.
SpeedyWeather.RingGrids.Matrix!
— MethodMatrix!(
MGs::Tuple{AbstractArray{T, 2}, OctahedralClenshawGrid}...;
quadrant_rotation,
matrix_quadrant
) -> Union{Tuple, AbstractMatrix}
Like Matrix!(::AbstractMatrix, ::OctahedralClenshawGrid)
but for simultaneous processing of tuples ((M1, G1), (M2, G2), ...)
with matrices Mi
and grids Gi
. All matrices and grids have to be of the same size respectively.
SpeedyWeather.RingGrids._scale_lat!
— Method_scale_lat!(
grid::AbstractGridArray{T, N, ArrayType} where {N, ArrayType<:AbstractArray{T, N}},
v::AbstractVector
) -> AbstractGridArray
Generic latitude scaling applied to A
in-place with latitude-like vector v
.
SpeedyWeather.RingGrids.anvil_average
— Methodanvil_average(a, b, c, d, Δab, Δcd, Δy) -> Any
The bilinear average of a, b, c, d which are values at grid points in an anvil-shaped configuration at location x, which is denoted by Δab, Δcd, Δy, the fraction of distances between a-b, c-d, and ab-cd, respectively. Note that a, c and b, d do not necessarily share the same longitude/x-coordinate. See schematic:
0..............1 # fraction of distance Δab between a, b
|< Δab >|
0^ a -------- o - b # anvil-shaped average of a, b, c, d at location x
.Δy |
. |
.v x
. |
1 c ------ o ---- d
|< Δcd >|
0...............1 # fraction of distance Δcd between c, d
^ fraction of distance Δy between a-b and c-d.
SpeedyWeather.RingGrids.average_on_poles
— Methodaverage_on_poles(
A::AbstractGridArray{NF<:Integer, 1, Array{NF<:Integer, 1}},
rings::Vector{<:UnitRange{<:Integer}}
) -> Tuple{Any, Any}
Method for A::Abstract{T<:Integer}
which rounds the averaged values to return the same number format NF
.
SpeedyWeather.RingGrids.average_on_poles
— Methodaverage_on_poles(
A::AbstractArray{NF<:AbstractFloat, 1},
rings::Vector{<:UnitRange{<:Integer}}
) -> Tuple{Any, Any}
Computes the average at the North and South pole from a given grid A
and it's precomputed ring indices rings
. The North pole average is an equally weighted average of all grid points on the northern-most ring. Similar for the South pole.
SpeedyWeather.RingGrids.check_inputs
— Methodcheck_inputs(data, nlat_half, rings, Grid) -> Any
True for data
, nlat_half
and rings
that all match in size to construct a grid of type Grid
.
SpeedyWeather.RingGrids.clenshaw_curtis_weights
— Methodclenshaw_curtis_weights(nlat_half::Integer) -> Any
The Clenshaw-Curtis weights for a Clenshaw grid (full or octahedral) of size nlathalf. Clenshaw-Curtis weights are of length nlat
, i.e. a vector for every latitude ring, pole to pole. `sum(clenshawcurtisweights(nlathalf))is always
2` as int0^π sin(x) dx = 2 (colatitudes), or equivalently int-pi/2^pi/2 cos(x) dx (latitudes).
Integration (and therefore the spectral transform) is exact (only rounding errors) when using Clenshaw grids provided that nlat >= 2(T + 1), meaning that a grid resolution of at least 128x64 (nlon x nlat) is sufficient for an exact transform with a T=31 spectral truncation.
SpeedyWeather.RingGrids.each_index_in_ring!
— Methodeach_index_in_ring!(
rings::Vector{<:UnitRange{<:Integer}},
Grid::Type{<:OctahedralGaussianArray},
nlat_half::Integer
)
precompute a Vector{UnitRange{Int}} to index grid points on every ring
j(elements of the vector) of
Gridat resolution
nlat_half. See
eachringand
eachgrid` for efficient looping over grid points.
SpeedyWeather.RingGrids.each_index_in_ring!
— Methodeach_index_in_ring!(
rings::Vector{<:UnitRange{<:Integer}},
Grid::Type{<:OctaminimalGaussianArray},
nlat_half::Integer
)
precompute a Vector{UnitRange{Int}} to index grid points on every ring
j(elements of the vector) of
Gridat resolution
nlat_half. See
eachringand
eachgrid` for efficient looping over grid points.
SpeedyWeather.RingGrids.each_index_in_ring
— Methodeach_index_in_ring(
Grid::Type{<:SpeedyWeather.RingGrids.AbstractFullGridArray},
j::Integer,
nlat_half::Integer
) -> Any
UnitRange
for every grid point of grid Grid
of resolution nlat_half
on ring j
(j=1
is closest ring around north pole, j=nlat
around south pole).
SpeedyWeather.RingGrids.each_index_in_ring
— Methodeach_index_in_ring(
grid::AbstractGridArray,
j::Integer
) -> Any
UnitRange to access data on grid grid
on ring j
.
SpeedyWeather.RingGrids.eachgrid
— Methodeachgrid(grid::AbstractGridArray) -> Any
CartesianIndices for the 2nd to last dimension of an AbstractGridArray, to be used like
for k in eachgrid(grid) for ring in eachring(grid) for ij in ring grid[ij, k]
SpeedyWeather.RingGrids.eachgridpoint
— Methodeachgridpoint(grid::AbstractGridArray) -> Base.OneTo
UnitRange to access each horizontal grid point on grid grid
. For a NxM
(N
horizontal grid points, M
vertical layers) OneTo(N)
is returned.
SpeedyWeather.RingGrids.eachgridpoint
— Methodeachgridpoint(
grid1::AbstractGridArray,
grids::Grid<:AbstractGridArray...
) -> Base.OneTo
Like eachgridpoint(::AbstractGridArray)
but checks for equal size between input arguments first.
SpeedyWeather.RingGrids.eachring
— Methodeachring(
grid1::AbstractGridArray,
grids::AbstractGridArray...
) -> Any
Same as eachring(grid)
but performs a bounds check to assess that all grids
according to grids_match
(non-parametric grid type, nlat_half and length).
SpeedyWeather.RingGrids.eachring
— Methodeachring(grid::AbstractGridArray) -> Any
Vector{UnitRange} rings
to loop over every ring of grid grid
and then each grid point per ring. To be used like
rings = eachring(grid)
for ring in rings
for ij in ring
grid[ij]
Accesses precomputed grid.rings
.
SpeedyWeather.RingGrids.eachring
— Methodeachring(
Grid::Type{<:AbstractGridArray},
nlat_half::Integer
) -> Any
Computes the ring indices i0:i1
for start and end of every longitudinal point on a given ring j
of Grid
at resolution nlat_half
. Used to loop over rings of a grid. These indices are also precomputed in every grid.rings
.
SpeedyWeather.RingGrids.equal_area_weights
— Methodequal_area_weights(
Grid::Type{<:AbstractGridArray},
nlat_half::Integer
) -> Any
The equal-area weights used for the HEALPix grids (original or OctaHEALPix) of size nlathalf. The weights are of length nlat
, i.e. a vector for every latitude ring, pole to pole. `sum(equalareaweights(nlathalf))is always
2` as int0^π sin(x) dx = 2 (colatitudes), or equivalently int-pi/2^pi/2 cos(x) dx (latitudes). Integration (and therefore the spectral transform) is not exact with these grids but errors reduce for higher resolution.
SpeedyWeather.RingGrids.extrema_in
— Methodextrema_in(v::AbstractVector, a::Real, b::Real) -> Any
For every element vᵢ in v does a<=vi<=b hold?
SpeedyWeather.RingGrids.full_array_type
— Methodfull_array_type(grid::AbstractGridArray) -> Any
Full grid array type for grid
. Always returns the N-dimensional *Array
not the two-dimensional (N=1
) *Grid
. For reduced grids the corresponding full grid that share the same latitudes.
SpeedyWeather.RingGrids.full_grid_type
— Methodfull_grid_type(grid::AbstractGridArray) -> Any
Full (horizontal) grid type for grid
. Always returns the two-dimensional (N=1
) *Grid
type. For reduced grids the corresponding full grid that share the same latitudes.
SpeedyWeather.RingGrids.gaussian_weights
— Methodgaussian_weights(nlat_half::Integer) -> Any
The Gaussian weights for a Gaussian grid (full or octahedral) of size nlathalf. Gaussian weights are of length nlat
, i.e. a vector for every latitude ring, pole to pole. `sum(gaussianweights(nlathalf))is always
2` as int0^π sin(x) dx = 2 (colatitudes), or equivalently int_-pi/2^pi/2 cos(x) dx (latitudes).
Integration (and therefore the spectral transform) is exact (only rounding errors) when using Gaussian grids provided that nlat >= 3(T + 1)/2, meaning that a grid resolution of at least 96x48 (nlon x nlat) is sufficient for an exact transform with a T=31 spectral truncation.
SpeedyWeather.RingGrids.get_colat
— Methodget_colat(grid::AbstractGridArray) -> Any
Colatitudes (radians) for meridional column (full grids only).
SpeedyWeather.RingGrids.get_colatlons
— Methodget_colatlons(grid::AbstractGridArray) -> Tuple{Any, Any}
Latitudes (in radians, 0-π) and longitudes (0 - 2π) for every (horizontal) grid point in grid
.
SpeedyWeather.RingGrids.get_lat
— Methodget_lat(grid::AbstractGridArray) -> Any
Latitude (radians) for each ring in grid
, north to south.
SpeedyWeather.RingGrids.get_latd
— Methodget_latd(grid::AbstractGridArray) -> Any
Latitude (degrees) for each ring in grid
, north to south.
SpeedyWeather.RingGrids.get_latdlonds
— Methodget_latdlonds(grid::AbstractGridArray) -> Tuple{Any, Any}
Latitudes (in degrees, -90˚-90˚N) and longitudes (0-360˚E) for every (horizontal) grid point in grid
. Ordered 0-360˚E then north to south.
SpeedyWeather.RingGrids.get_latlons
— Methodget_latlons(grid::AbstractGridArray) -> Tuple{Any, Any}
Latitudes (in radians, 0-2π) and longitudes (-π/2 - π/2) for every (horizontal) grid point in grid
.
SpeedyWeather.RingGrids.get_lon
— Methodget_lon(grid::AbstractGridArray) -> Any
Longitude (radians) for meridional column (full grids only).
SpeedyWeather.RingGrids.get_lond
— Methodget_lond(grid::AbstractGridArray) -> Any
Longitude (degrees) for meridional column (full grids only).
SpeedyWeather.RingGrids.get_nlat
— Methodget_nlat(grid::AbstractGridArray) -> Any
Get number of latitude rings, pole to pole.
SpeedyWeather.RingGrids.get_nlat_half
— Methodget_nlat_half(grid::AbstractGridArray) -> Any
Resolution paraemeters nlat_half
of a grid
. Number of latitude rings on one hemisphere, Equator included.
SpeedyWeather.RingGrids.get_nlon_per_ring
— Methodget_nlon_per_ring(
Grid::Type{<:HEALPixArray},
nlat_half::Integer,
j::Integer
) -> Any
Number of longitude points for ring j
on Grid
of resolution nlat_half
.
SpeedyWeather.RingGrids.get_nlons
— Methodget_nlons(
Grid::Type{<:AbstractGridArray},
nlat_half::Integer;
both_hemispheres
) -> Any
Returns a vector nlons
for the number of longitude points per latitude ring, north to south. Provide grid Grid
and its resolution parameter nlat_half
. For keyword argument both_hemispheres=false
only the northern hemisphere (incl Equator) is returned.
SpeedyWeather.RingGrids.get_npoints
— Methodget_npoints(grid::AbstractGridArray) -> Any
Total number of grid points in all dimensions of grid
. Equivalent to length of the underlying data array.
SpeedyWeather.RingGrids.get_npoints2D
— Methodget_npoints2D(grid::AbstractGridArray) -> Any
Number of grid points in the horizontal dimension only, even if grid
is N-dimensional.
SpeedyWeather.RingGrids.grid_cell_average!
— Methodgrid_cell_average!(
output::AbstractGrid,
input::SpeedyWeather.RingGrids.AbstractFullGrid
) -> AbstractGrid
Averages all grid points in input
that are within one grid cell of output
with coslat-weighting. The output grid cell boundaries are assumed to be rectangles spanning half way to adjacent longitude and latitude points.
SpeedyWeather.RingGrids.grid_cell_average
— Methodgrid_cell_average(
Grid::Type{<:AbstractGrid},
nlat_half::Integer,
input::SpeedyWeather.RingGrids.AbstractFullGrid
) -> Any
Averages all grid points in input
that are within one grid cell of output
with coslat-weighting. The output grid cell boundaries are assumed to be rectangles spanning half way to adjacent longitude and latitude points.
SpeedyWeather.RingGrids.grids_match
— Methodgrids_match(
A::AbstractGridArray,
B::AbstractGridArray;
horizontal_only,
vertical_only
) -> Any
True if both A
and B
are of the same nonparametric grid type (e.g. OctahedralGaussianArray, regardless type parameter T
or underyling array type ArrayType
) and of same resolution (nlat_half
) and total grid points (length
). Sizes of (4,)
and (4,1)
would match for example, but (8,1)
and (4,2)
would not (nlat_half
not identical).
SpeedyWeather.RingGrids.grids_match
— Methodgrids_match(
A::AbstractGridArray,
B::AbstractGridArray...;
kwargs...
) -> Any
True if all grids A, B, C, ...
provided as arguments match according to grids_match
wrt to A
(and therefore all).
SpeedyWeather.RingGrids.horizontal_grid_type
— Methodhorizontal_grid_type(grid::AbstractGridArray) -> Any
The two-dimensional (N=1
) *Grid
for grid
, which can be an N-dimensional *GridArray
.
SpeedyWeather.RingGrids.isdecreasing
— Methodisdecreasing(x::AbstractVector) -> Bool
Check whether elements of a vector v
are strictly decreasing.
SpeedyWeather.RingGrids.isincreasing
— Methodisincreasing(x::AbstractVector) -> Bool
Check whether elements of a vector v
are strictly increasing.
SpeedyWeather.RingGrids.matrix_size
— Methodmatrix_size(grid::AbstractGridArray) -> Tuple{Int64, Int64}
Size of the matrix of the horizontal grid if representable as such (not all grids).
SpeedyWeather.RingGrids.nlat_odd
— Methodnlat_odd(grid::AbstractGridArray) -> Any
True for a grid
that has an odd number of latitude rings nlat
(both hemispheres).
SpeedyWeather.RingGrids.nonparametric_type
— Methodnonparametric_type(grid::AbstractGridArray) -> Any
For any instance of AbstractGridArray
type its n-dimensional type (*Grid{T, N, ...} returns *Array) but without any parameters {T, N, ArrayType}
SpeedyWeather.RingGrids.npoints_added_per_ring
— Methodnpoints_added_per_ring(
_::Type{<:OctahedralGaussianArray}
) -> Int64
[EVEN MORE EXPERIMENTAL] number of longitude points added (removed) for every ring towards the Equator (on the southern hemisphere towards the south pole).
SpeedyWeather.RingGrids.npoints_added_per_ring
— Methodnpoints_added_per_ring(
_::Type{<:OctaminimalGaussianArray}
) -> Int64
[EVEN MORE EXPERIMENTAL] number of longitude points added (removed) for every ring towards the Equator (on the southern hemisphere towards the south pole).
SpeedyWeather.RingGrids.npoints_pole
— Methodnpoints_pole(_::Type{<:OctahedralGaussianArray}) -> Int64
[EXPERIMENTAL] additional number of longitude points on the first and last ring. Change to 0 to start with 4 points on the first ring.
SpeedyWeather.RingGrids.npoints_pole
— Methodnpoints_pole(_::Type{<:OctaminimalGaussianArray}) -> Int64
[EXPERIMENTAL] additional number of longitude points on the first and last ring. Change to 0 to start with 4 points on the first ring.
SpeedyWeather.RingGrids.nside_healpix
— Methodnside_healpix(nlat_half::Integer) -> Any
The original Nside
resolution parameter of the HEALPix grids. The number of grid points on one side of each (square) face. While we use nlat_half
across all ring grids, this function translates this to Nside. Even nlat_half
only.
SpeedyWeather.RingGrids.rotate_matrix_indices_180
— Methodrotate_matrix_indices_180(
i::Integer,
j::Integer,
s::Integer
) -> Tuple{Any, Any}
Rotate indices i, j
of a square matrix of size s x s by 180˚.
SpeedyWeather.RingGrids.rotate_matrix_indices_270
— Methodrotate_matrix_indices_270(
i::Integer,
j::Integer,
s::Integer
) -> Tuple{Integer, Any}
Rotate indices i, j
of a square matrix of size s x s anti-clockwise by 270˚.
SpeedyWeather.RingGrids.rotate_matrix_indices_90
— Methodrotate_matrix_indices_90(
i::Integer,
j::Integer,
s::Integer
) -> Tuple{Any, Integer}
Rotate indices i, j
of a square matrix of size s x s anti-clockwise by 90˚.
SpeedyWeather.RingGrids.spherical_distance
— Methodspherical_distance(
Formula::Type{<:SpeedyWeather.RingGrids.AbstractSphericalDistance},
args...;
kwargs...
) -> Any
Spherical distance, or great-circle distance, between two points lonlat1
and lonlat2
using the Formula
(default Haversine
).
SpeedyWeather.RingGrids.whichring
— Methodwhichring(
ij::Integer,
rings::Vector{UnitRange{Int64}}
) -> Int64
Obtain ring index j
from gridpoint ij
and rings
describing rind indices as obtained from eachring(::Grid)