Rasters.jl
Rasters
— ModuleRasters
Rasters.jl defines common types and methods for reading, writing and manipulating rasterized spatial data.
These currently include raster arrays like GeoTIFF and NetCDF, R grd files, multi-layered stacks, and multi-file series of arrays and stacks.
A RasterStack of EarthEnv HabitatHeterogeneity layers, trimmed to Australia and plotted with Plots.jl
:warning: Packages extensions and Rasters 0.8 and onwards
On Julia 1.9 we can put additional packages in extensions, so the code only loads when you load a specific package. Rasters.jl was always intended to work like this, and its finally possible. This reduced package using
time from many seconds to well under a second.
But, it means you have to manually load packages you need for each backend or additional functionality.
For example, to use the GDAL backend, and download files, you now need to do:
using Rasters, ArchGDAL, RasterDataSources
where previously it was just using Rasters
.
Sources and packages needed:
:gdal
:using ArchGDAL
:netcdf
:using NCDatasets
:grd
: built-in.:smap
:using HDF5
:grib
: not yet finished.
Other functionality in extensions:
- Raster data downloads, like
Worldclim{Climate}
:using RasterDataSources
- Makie plots:
using Makie
- Coordinate transformations for gdal rasters:
using CoordinateTransformations
Quick start
Install the package by typing:
]
+add Rasters
using Rasters
Using Rasters
to read GeoTiff or NetCDF files will output something similar to the following toy examples. This is possible because Rasters.jl extends DimensionalData.jl so that spatial data can be indexed using named dimensions like X
, Y
and Ti
(time) and e.g. spatial coordinates.
using Rasters, Dates
+lon, lat = X(25:1:30), Y(25:1:30)
+ti = Ti(DateTime(2001):Month(1):DateTime(2002))
+ras = Raster(rand(lon, lat, ti)) # this generates random numbers with the dimensions given
6×6×13 Raster{Float64,3} with dimensions:
+ X Sampled{Int64} 25:1:30 ForwardOrdered Regular Points,
+ Y Sampled{Int64} 25:1:30 ForwardOrdered Regular Points,
+ Ti Sampled{DateTime} DateTime("2001-01-01T00:00:00"):Month(1):DateTime("2002-01-01T00:00:00") ForwardOrdered Regular Points
+extent: Extent(X = (25, 30), Y = (25, 30), Ti = (DateTime("2001-01-01T00:00:00"), DateTime("2002-01-01T00:00:00")))
+missingval: missing
+values: [:, :, 1]
+ 25 26 27 28 29 30
+ 25 0.9063 0.427328 0.0320967 0.297023 0.0571002 0.891377
+ 26 0.443494 0.867547 0.350546 0.150155 0.24565 0.711039
+ 27 0.745673 0.0991336 0.930332 0.893537 0.805931 0.360583
+ 28 0.512083 0.125287 0.959434 0.354868 0.337824 0.259563
+ 29 0.253849 0.692209 0.774092 0.131798 0.823656 0.390013
+ 30 0.334152 0.136551 0.183555 0.941133 0.450484 0.461862
+[and 12 more slices...]
Getting the lookup
array from dimensions
lon = lookup(ras, X) # if X is longitude
+lat = lookup(ras, Y) # if Y is latitude
Sampled{Int64} ForwardOrdered Regular Points
+wrapping: 25:1:30
Select by index
Selecting a time slice by index
is done via
ras[Ti(1)]
6×6 Raster{Float64,2} with dimensions:
+ X Sampled{Int64} 25:1:30 ForwardOrdered Regular Points,
+ Y Sampled{Int64} 25:1:30 ForwardOrdered Regular Points
+and reference dimensions:
+ Ti Sampled{DateTime} DateTime("2001-01-01T00:00:00"):Month(1):DateTime("2001-01-01T00:00:00") ForwardOrdered Regular Points
+extent: Extent(X = (25, 30), Y = (25, 30))
+missingval: missing
+values: 25 26 27 28 29 30
+ 25 0.9063 0.427328 0.0320967 0.297023 0.0571002 0.891377
+ 26 0.443494 0.867547 0.350546 0.150155 0.24565 0.711039
+ 27 0.745673 0.0991336 0.930332 0.893537 0.805931 0.360583
+ 28 0.512083 0.125287 0.959434 0.354868 0.337824 0.259563
+ 29 0.253849 0.692209 0.774092 0.131798 0.823656 0.390013
+ 30 0.334152 0.136551 0.183555 0.941133 0.450484 0.461862
ras[Ti=1]
6×6 Raster{Float64,2} with dimensions:
+ X Sampled{Int64} 25:1:30 ForwardOrdered Regular Points,
+ Y Sampled{Int64} 25:1:30 ForwardOrdered Regular Points
+and reference dimensions:
+ Ti Sampled{DateTime} DateTime("2001-01-01T00:00:00"):Month(1):DateTime("2001-01-01T00:00:00") ForwardOrdered Regular Points
+extent: Extent(X = (25, 30), Y = (25, 30))
+missingval: missing
+values: 25 26 27 28 29 30
+ 25 0.9063 0.427328 0.0320967 0.297023 0.0571002 0.891377
+ 26 0.443494 0.867547 0.350546 0.150155 0.24565 0.711039
+ 27 0.745673 0.0991336 0.930332 0.893537 0.805931 0.360583
+ 28 0.512083 0.125287 0.959434 0.354868 0.337824 0.259563
+ 29 0.253849 0.692209 0.774092 0.131798 0.823656 0.390013
+ 30 0.334152 0.136551 0.183555 0.941133 0.450484 0.461862
or and interval of indices using the syntax =a:b
or (a:b)
ras[Ti(1:10)]
6×6×10 Raster{Float64,3} with dimensions:
+ X Sampled{Int64} 25:1:30 ForwardOrdered Regular Points,
+ Y Sampled{Int64} 25:1:30 ForwardOrdered Regular Points,
+ Ti Sampled{DateTime} DateTime("2001-01-01T00:00:00"):Month(1):DateTime("2001-10-01T00:00:00") ForwardOrdered Regular Points
+extent: Extent(X = (25, 30), Y = (25, 30), Ti = (DateTime("2001-01-01T00:00:00"), DateTime("2001-10-01T00:00:00")))
+missingval: missing
+values: [:, :, 1]
+ 25 26 27 28 29 30
+ 25 0.9063 0.427328 0.0320967 0.297023 0.0571002 0.891377
+ 26 0.443494 0.867547 0.350546 0.150155 0.24565 0.711039
+ 27 0.745673 0.0991336 0.930332 0.893537 0.805931 0.360583
+ 28 0.512083 0.125287 0.959434 0.354868 0.337824 0.259563
+ 29 0.253849 0.692209 0.774092 0.131798 0.823656 0.390013
+ 30 0.334152 0.136551 0.183555 0.941133 0.450484 0.461862
+[and 9 more slices...]
Select by value
ras[Ti=At(DateTime(2001))]
6×6 Raster{Float64,2} with dimensions:
+ X Sampled{Int64} 25:1:30 ForwardOrdered Regular Points,
+ Y Sampled{Int64} 25:1:30 ForwardOrdered Regular Points
+and reference dimensions:
+ Ti Sampled{DateTime} DateTime("2001-01-01T00:00:00"):Month(1):DateTime("2001-01-01T00:00:00") ForwardOrdered Regular Points
+extent: Extent(X = (25, 30), Y = (25, 30))
+missingval: missing
+values: 25 26 27 28 29 30
+ 25 0.9063 0.427328 0.0320967 0.297023 0.0571002 0.891377
+ 26 0.443494 0.867547 0.350546 0.150155 0.24565 0.711039
+ 27 0.745673 0.0991336 0.930332 0.893537 0.805931 0.360583
+ 28 0.512083 0.125287 0.959434 0.354868 0.337824 0.259563
+ 29 0.253849 0.692209 0.774092 0.131798 0.823656 0.390013
+ 30 0.334152 0.136551 0.183555 0.941133 0.450484 0.461862
More options are available, like Near
, Contains
and Where
. For more details go here.
Dimensions can also be used in most Base
and Statistics
methods like mean
and reduce
where dims
arguments are required. Much of the behaviour is covered in the DimensionalData docs.
See the docs for more details and examples for Rasters.jl.
Data-source abstraction
Rasters provides a standardised interface that allows many source data types to be used with identical syntax.
- Scripts and packages building on Rasters.jl can treat
Raster
,RasterStack
, andRasterSeries
as black boxes.- The data could hold GeoTiff or NetCDF files,
Array
s in memory orCuArray
s on the GPU - they will all behave in the same way. RasterStack
can be backed by a Netcdf or HDF5 file, or aNamedTuple
ofRaster
holding.tif
files, or allRaster
in memory.- Users do not have to deal with the specifics of spatial file types.
- The data could hold GeoTiff or NetCDF files,
Projected
lookups with Cylindrical projections can by indexed using other Cylindrical projections by setting themappedcrs
keyword on construction. You don't need to know the underlying projection, the conversion is handled automatically. This means lat/lonEPSG(4326)
can be used seamlessly if you need that.
Bugs, errors and making issues for Rasters.jl
Raster data is complicated and there are many places for subtle or not-so-subtle bugs to creep in.
We need bug reports to reduce how often they occur over time. But also, we need issues that are easy to reproduce or it isn't practically possible to fix them.
Because there are so many raster file types and variations of them, most of the time we need the exact file that caused your problem to know how to fix it, and be sure that we have actually fixed it when we are done. So fixing a Rasters.jl bug nearly always involves downloading some file and running some code that breaks with it (if you can trigger the bug without a file, thats great! but its not always possible).
To make an issue we can fix quickly (or at all) there are three key steps:
- Include the file in an accessible place on web without autentication or any other work on our part, so we can just get it and find your bug. You can put it on a file hosting platform (e.g. google drive, drop box, whatever you use) and share the url.
- Add a minimum working example to the issue template that first downloads the file, then runs the function that triggers the bug.
- Paste the complete stack trace of the error it produces, right to the bottom, into the issue template. Then we can be sure we reproduced the same problem.
Good issues are really appreciated, but they do take just a little extra effort with Rasters.jl because of this need for files.
Common Applications
Subsetting an object
Regular getindex
(e.g. A[1:100, :]
) and view
work on all objects just as with an Array
. view
is always lazy, and reads from disk are deferred until getindex
is used. DimensionalData.jl Dimension
s and Selector
s are the other way to subset an object, making use of the objects index to find values at e.g. certain X/Y coordinates. The available selectors are listed here:
At(x) | get the index exactly matching the passed in value(s). |
Near(x) | get the closest index to the passed in value(s). |
Where(f::Function) | filter the array axis by a function of the dimension index values. |
a..b /Between(a, b) | get all indices between two values, excluding the high value. |
Contains(x) | get indices where the value x falls within an interval. |
Use the ..
selector to take a view
of madagascar:
using Rasters, RasterDataSources, ArchGDAL, Plots
+A = Raster(WorldClim{BioClim}, 5)
+madagascar = view(A, X(43.25 .. 50.48), Y(-25.61 .. -12.04)) # Note the space between .. -12
+plot(madagascar)
+
+Methods that change the reslolution or extent of an object
Click through to the function documentation for more in-depth descriptions and examples.
aggregate | aggregate data by the same or different amounts for each axis. |
disaggregate | similarly disaggregate data. |
mosaic | join rasters covering different extents into a single array or file. |
crop | shrink objects to specific dimension sizes or the extent of another object. |
extend | extend objects to specific dimension sizes or the extent of another object. |
trim | trims areas of missing values for arrays and across stack layers. |
resample | resample data to a different size and projection, or snap to another object. |
warp | use gdalwarp on any object, e.g. a multidimensional NetCDF stack. |
Methods that change an objects values:
Note that most regular Julia methods, such as replace
, work as for a standard Array
. These additional methods are commonly required in GIS applications.
classify | classify values into categories. |
mask | mask an object by a polygon or Raster along X/Y , or other dimensions. |
replace_missing | replace all missing values in an object and update missingval . |
Point, polygon and table operation
rasterize | rasterize points and geometries. |
coverage | get the fraction of each pixel covered by geometries. |
extract | extract values from points or geometries. |
zonal | calculate zonal statistics for an object masked by geometries. |
Methods to load, write and modify data sources:
modify | replace the data in objects. Useful to e.g. move objects to/from a GPU. |
read | read data to memory if it is on disk. |
read! | read data to predefined memory. |
open | open the underlying data for manually reading or writing. |
write | write objects to file. |
Altering and summarising arrays and stacks with regular julia methods
Most base methods work as for regular julia Array
s, such as reverse
and rotations like rotl90
. Base, statistics and linear algebra methods like mean
that take a dims
argument can also use the dimension name. To take the mean over the time dimension:
mean(A, dims=Ti)
broadcast
works lazily from disk when lazy=true
, and is only applied when data is directly indexed. Adding a dot to any function will use broadcast over a Raster
just like an Array
.
Broadcasting
For a disk-based array A
, this will only be applied when indexing occurs or when we read
the array.
A .*= 2
To broadcast directly to disk, we need to open the file in write mode:
open(Raster(filename); write=true) do O
+ O .*= 2
+end
To broadcast over a RasterStack
use map
, which applies a function to the raster layers of the stack.
newstack = map(stack) do raster
+ raster .* 2
+end
Modifying object properties
rebuild
can be used to modify the fields of an object, generating a new object (but possibly holding the same arrays or files).
If you know that a file had an incorrectly specified missing value, you can do:
rebuild(A; missingval=-9999)
(replace_missing
will actualy replace the current values)
Or if you need to change the name of the layer:
rebuild(A; name=:temperature)
set
can be used to modify the nested properties of an objects dimensions, that are more difficult to change with rebuild
. set
works on the principal that dimension properties can only be in one specific field, so we generally don't have to specify which one it is. set
will also try to update anything affected by a change you make.
This will set the X
axis to specify points, instead of intervals:
set(A, X => Points)
We can also reassign dimensions, here X
becomes Z
:
set(A, X => Z)
setcrs(A, crs)
and setmappedcrs(A, crs)
will set the crs value/s of an object to any GeoFormat
from GeoFormatTypes.jl.
Examples and Plotting
Plots.jl and Makie.jl are fully supported by Rasters.jl, with recipes for plotting Raster
and RasterStack
provided. plot
will plot a heatmap with axes matching dimension values. If mappedcrs
is used, converted values will be shown on axes instead of the underlying crs
values. contourf
will similarly plot a filled contour plot.
Pixel resolution is limited to allow loading very large files quickly. max_res
specifies the maximum pixel resolution to show on the longest axis of the array. It can be set manually to change the resolution (e.g. for large or high-quality plots):
using Rasters, RasterDataSources, ArchGDAL, Plots
+A = Raster(WorldClim{BioClim}, 5)
+plot(A; max_res=3000)
For Makie, plot
functions in a similar way. plot
will only accept two-dimensional rasters. You can invoke contour
, contourf
, heatmap
, surface
or any Makie plotting function which supports surface-like data on a 2D raster.
To obtain tiled plots for 3D rasters and RasterStacks, use the function Rasters.rplot([gridposition], raster; kw_args...)
. This is an unexported function, since we're not sure how the API will change going forward.
using Rasters, CairoMakie, RasterDataSources, ArchGDAL
+A = Raster(WorldClim{BioClim}, 5)
+Makie.plot(A)
Loading and plotting data
Our first example simply loads a file from disk and plots it.
This netcdf file only has one layer, if it has more we could use RasterStack
instead.
using Rasters, NCDatasets, Plots
+url = "https://www.unidata.ucar.edu/software/netcdf/examples/tos_O1_2001-2002.nc";
+filename = download(url, "tos_O1_2001-2002.nc");
+A = Raster(filename)
180×170×24 Raster{Union{Missing, Float32},3} tos with dimensions:
+ X Mapped{Float64} Float64[1.0, 3.0, …, 357.0, 359.0] ForwardOrdered Explicit Intervals crs: EPSG mappedcrs: EPSG,
+ Y Mapped{Float64} Float64[-79.5, -78.5, …, 88.5, 89.5] ForwardOrdered Explicit Intervals crs: EPSG mappedcrs: EPSG,
+ Ti Sampled{DateTime360Day} DateTime360Day[DateTime360Day(2001-01-16T00:00:00), …, DateTime360Day(2002-12-16T00:00:00)] ForwardOrdered Explicit Intervals
+extent: Extent(X = (-0.0, 360.0), Y = (-80.0, 90.0), Ti = (DateTime360Day(2001-01-01T00:00:00), DateTime360Day(2003-01-01T00:00:00)))missingval: missingcrs: EPSG:4326
+mappedcrs: EPSG:4326
+parent:
+[:, :, 1]
+ -79.5 -78.5 … 86.5 87.5 88.5 89.5
+ 1.0 missing missing 271.43 271.437 271.445 271.459
+ 3.0 missing missing 271.431 271.438 271.445 271.459
+ 5.0 missing missing 271.431 271.438 271.445 271.459
+ 7.0 missing missing 271.431 271.439 271.446 271.459
+ ⋮ ⋱ ⋮
+ 351.0 missing missing 271.43 271.435 271.445 271.459
+ 353.0 missing missing 271.43 271.436 271.445 271.459
+ 355.0 missing missing 271.43 271.436 271.445 271.459
+ 357.0 missing missing 271.43 271.437 271.445 271.459
+ 359.0 missing missing … 271.43 271.437 271.445 271.459
+[and 23 more slices...]
Objects with Dimensions other than X
and Y
will produce multi-pane plots. Here we plot every third month in the first year in one plot:
A[Ti=1:3:12] |> plot
+
+Now plot the ocean temperatures around the Americas in the first month of 2001. Notice we are using lat/lon coordinates and date/time instead of regular indices. The time dimension uses DateTime360Day
, so we need to load CFTime.jl to index it with Near
.
using CFTime
+A[Ti(Near(DateTime360Day(2001, 01, 17))), Y(-60.0 .. 90.0), X(45.0 .. 190.0)] |> plot
+
+Now get the mean over the timespan, then save it to disk, and plot it as a filled contour:
Other plot functions and sliced objects that have only one X
/Y
/Z
dimension fall back to generic DimensionalData.jl plotting, which will still correctly label plot axes.
using Statistics
+# Take the mean
+mean_tos = mean(A; dims=Ti)
180×170×1 Raster{Union{Missing, Float32},3} tos with dimensions:
+ X Mapped{Float64} Float64[1.0, 3.0, …, 357.0, 359.0] ForwardOrdered Explicit Intervals crs: EPSG mappedcrs: EPSG,
+ Y Mapped{Float64} Float64[-79.5, -78.5, …, 88.5, 89.5] ForwardOrdered Explicit Intervals crs: EPSG mappedcrs: EPSG,
+ Ti Sampled{DateTime360Day} DateTime360Day(2002-01-16T00:00:00):Millisecond(2592000000):DateTime360Day(2002-01-16T00:00:00) ForwardOrdered Explicit Intervals
+extent: Extent(X = (-0.0, 360.0), Y = (-80.0, 90.0), Ti = (DateTime360Day(2001-01-01T00:00:00), DateTime360Day(2003-01-01T00:00:00)))missingval: missingcrs: EPSG:4326
+mappedcrs: EPSG:4326
+parent:
+[:, :, 1]
+ -79.5 -78.5 … 86.5 87.5 88.5 89.5
+ 1.0 missing missing 271.427 271.434 271.443 271.454
+ 3.0 missing missing 271.427 271.434 271.443 271.454
+ 5.0 missing missing 271.427 271.434 271.443 271.454
+ 7.0 missing missing 271.427 271.435 271.444 271.454
+ ⋮ ⋱ ⋮
+ 351.0 missing missing 271.427 271.433 271.443 271.454
+ 353.0 missing missing 271.427 271.433 271.443 271.454
+ 355.0 missing missing 271.427 271.433 271.443 271.454
+ 357.0 missing missing 271.427 271.433 271.443 271.454
+ 359.0 missing missing … 271.427 271.433 271.443 271.454
Plot a contour plot
contourf(mean_tos; dpi=300, size=(800, 400))
+
+Write the mean values to disk
write("mean_tos.nc", mean_tos)
"mean_tos.nc"
Plotting recipes in DimensionalData.jl are the fallback for Rasters.jl when the object doesn't have 2 X
/Y
/Z
dimensions, or a non-spatial plot command is used. So (as a random example) we could plot a transect of ocean surface temperature at 20 degree latitude :
A[Y(Near(20.0)), Ti(1)] |> plot
+
+A basic species distribution modelling workflow
Load occurrences for the Mountain Pygmy Possum using GBIF.jl
using Rasters, RasterDataSources, ArchGDAL, GBIF2, Plots
+records = GBIF2.occurrence_search("Burramys parvus"; limit=300)
300-element GBIF2.Table{GBIF2.Occurrence, JSON3.Array{JSON3.Object, Vector{UInt8}, SubArray{UInt64, 1, Vector{UInt64}, Tuple{UnitRange{Int64}}, true}}}┌──────────────────┬─────────────────┬─────────┬─────────┬─────────┬────────────
+│ decimalLongitude │ decimalLatitude │ year │ month │ day │ kingdom ⋯
+│ Float64? │ Float64? │ Int64? │ Int64? │ Int64? │ String? ⋯
+├──────────────────┼─────────────────┼─────────┼─────────┼─────────┼────────────
+│ missing │ missing │ missing │ missing │ missing │ Animalia ⋯
+│ missing │ missing │ missing │ missing │ missing │ Animalia ⋯
+│ missing │ missing │ 2021 │ 1 │ 6 │ Animalia ⋯
+│ missing │ missing │ missing │ missing │ missing │ Animalia ⋯
+│ missing │ missing │ missing │ missing │ missing │ Animalia ⋯
+│ missing │ missing │ missing │ missing │ missing │ Animalia ⋯
+│ 148.391 │ -36.3036 │ 2015 │ 11 │ 15 │ Animalia ⋯
+│ 148.333 │ -36.4333 │ 2011 │ 11 │ 21 │ Animalia ⋯
+│ 148.396 │ -36.3818 │ 2016 │ 11 │ 15 │ Animalia ⋯
+│ missing │ missing │ 1975 │ 3 │ 26 │ Animalia ⋯
+│ 147.096 │ -36.9357 │ 2020 │ 2 │ 10 │ Animalia ⋯
+│ 148.329 │ -36.4317 │ 2016 │ 1 │ 3 │ Animalia ⋯
+│ 148.236 │ -36.5249 │ 2012 │ 11 │ 23 │ Animalia ⋯
+│ missing │ missing │ missing │ missing │ missing │ Animalia ⋯
+│ missing │ missing │ missing │ missing │ missing │ Animalia ⋯
+│ ⋮ │ ⋮ │ ⋮ │ ⋮ │ ⋮ │ ⋮ ⋱
+└──────────────────┴─────────────────┴─────────┴─────────┴─────────┴────────────
+ 77 columns and 285 rows omitted
+
Extract the longitude/latitude value to a Vector
of points (a Tuple
counts as a (x, y)
point in GeoInterface.jl):
coords = [(r.decimalLongitude, r.decimalLatitude) for r in records if !ismissing(r.decimalLatitude)]
253-element Vector{Tuple{Float64, Float64}}:
+ (148.391097, -36.30362)
+ (148.332969, -36.433349)
+ (148.396453, -36.381847)
+ (147.096394, -36.935687)
+ (148.328896, -36.431684)
+ (148.235596, -36.524924)
+ (148.240881, -36.400058)
+ (148.347186, -36.504673)
+ (151.250866, -33.831883)
+ (147.1, -37.0)
+ ⋮
+ (148.248473, -36.424117)
+ (148.368407, -36.413914)
+ (148.256519, -36.442973)
+ (148.467196, -36.000859)
+ (148.40889, -36.006805)
+ (148.40889, -36.006805)
+ (148.40889, -36.006805)
+ (148.467196, -36.000859)
+ (148.392241, -35.951583)
Get BioClim layers and subset to south-east Australia
A = RasterStack(WorldClim{BioClim}, (1, 3, 7, 12))
+se_aus = A[X(138 .. 155), Y(-40 .. -25)]
RasterStack with dimensions:
+ X Projected{Float64} LinRange{Float64}(138.0, 154.833, 102) ForwardOrdered Regular Intervals crs: WellKnownText,
+ Y Projected{Float64} LinRange{Float64}(-25.1667, -39.8333, 89) ReverseOrdered Regular Intervals crs: WellKnownText
+and 4 layers:
+ :bio1 Float32 dims: X, Y (102×89)
+ :bio3 Float32 dims: X, Y (102×89)
+ :bio7 Float32 dims: X, Y (102×89)
+ :bio12 Float32 dims: X, Y (102×89)
+
Plot BioClim predictors and scatter occurrence points on all subplots
p = plot(se_aus);
+kw = (legend=:none, opacity=0.5, markershape=:cross, markercolor=:black)
+foreach(i -> scatter!(p, coords; subplot=i, kw...), 1:4)
+display(p)
Then extract predictor variables and write to CSV.
using CSV
+predictors = collect(extract(se_aus, coords))
+CSV.write("burramys_parvus_predictors.csv", predictors)
"burramys_parvus_predictors.csv"
Or convert them to a DataFrame
.
using DataFrames
+df = DataFrame(predictors)
+df[1:5, :]
Row | geometry | bio1 | bio3 | bio7 | bio12 |
---|---|---|---|---|---|
Tuple… | Float32 | Float32 | Float32 | Float32 | |
1 | (148.391, -36.3036) | 6.1707 | 41.1198 | 23.4645 | 1692.0 |
2 | (148.333, -36.4333) | 7.83572 | 41.5975 | 23.5028 | 1500.0 |
3 | (148.396, -36.3818) | 6.88158 | 42.2681 | 23.133 | 1544.0 |
4 | (147.096, -36.9357) | 9.40835 | 40.7905 | 23.0895 | 1292.0 |
5 | (148.329, -36.4317) | 7.83572 | 41.5975 | 23.5028 | 1500.0 |
Polygon masking, mosaic and plot
In this example we will mask
the Scandinavian countries with border polygons, then mosaic
together to make a single plot.
First, get the country boundary shape files using GADM.jl.
using Rasters, RasterDataSources, ArchGDAL, Shapefile, Plots, Dates, Downloads, NCDatasets
+
+# Download the shapefile
+shapefile_url = "https://github.com/nvkelso/natural-earth-vector/raw/master/10m_cultural/ne_10m_admin_0_countries.shp"
+shapefile_name = "boundary_lines.shp"
+Downloads.download(shapefile_url, shapefile_name)
+
+# Load using Shapefile.jl
+shapes = Shapefile.Handle(shapefile_name)
+denmark_border = shapes.shapes[71]
+norway_border = shapes.shapes[53]
+sweden_border = shapes.shapes[54]
Polygon(4665 Points)
Then load raster data. We load some worldclim layers using RasterDataSources
via Rasters.jl:
climate = RasterStack(WorldClim{Climate}, (:tmin, :tmax, :prec, :wind); month=July)
RasterStack with dimensions:
+ X Projected{Float64} LinRange{Float64}(-180.0, 179.833, 2160) ForwardOrdered Regular Intervals crs: WellKnownText,
+ Y Projected{Float64} LinRange{Float64}(89.8333, -90.0, 1080) ReverseOrdered Regular Intervals crs: WellKnownText
+and 4 layers:
+ :tmin Float32 dims: X, Y (2160×1080)
+ :tmax Float32 dims: X, Y (2160×1080)
+ :prec Int16 dims: X, Y (2160×1080)
+ :wind Float32 dims: X, Y (2160×1080)
+
mask
Denmark, Norway and Sweden from the global dataset using their border polygon, then trim the missing values. We pad trim
with a 10 pixel margin.
mask_trim(climate, poly) = trim(mask(climate; with=poly); pad=10)
+
+denmark = mask_trim(climate, denmark_border)
+norway = mask_trim(climate, norway_border)
+sweden = mask_trim(climate, sweden_border)
RasterStack with dimensions:
+ X Projected{Float64} LinRange{Float64}(9.5, 25.6667, 98) ForwardOrdered Regular Intervals crs: WellKnownText,
+ Y Projected{Float64} LinRange{Float64}(70.5, 53.6667, 102) ReverseOrdered Regular Intervals crs: WellKnownText
+and 4 layers:
+ :tmin Float32 dims: X, Y (98×102)
+ :tmax Float32 dims: X, Y (98×102)
+ :prec Int16 dims: X, Y (98×102)
+ :wind Float32 dims: X, Y (98×102)
+
Plotting
First define a function to add borders to all subplots.
function borders!(p, poly)
+ for i in 1:length(p)
+ plot!(p, poly; subplot=i, fillalpha=0, linewidth=0.6)
+ end
+ return p
+end
borders! (generic function with 1 method)
Now we can plot the individual countries.
dp = plot(denmark)
+borders!(dp, denmark_border)
+
+sp = plot(sweden)
+borders!(sp, sweden_border)
+
+np = plot(norway)
+borders!(np, norway_border)
+
+The Norway shape includes a lot of islands. Lets crop them out using ..
intervals:
norway_region = climate[X(0..40), Y(55..73)]
+plot(norway_region)
+
+And mask it with the border again:
norway = mask_trim(norway_region, norway_border)
+np = plot(norway)
+borders!(np, norway_border)
+
+Now we can combine the countries into a single raster using mosaic
. first
will take the first value if/when there is an overlap.
scandinavia = mosaic(first, denmark, norway, sweden)
RasterStack with dimensions:
+ X Projected{Float64} 3.1666666666666443:0.16666666666666666:32.49999999999998 ForwardOrdered Regular Intervals crs: WellKnownText,
+ Y Projected{Float64} 72.66666666666666:-0.16666666666666666:52.99999999999999 ReverseOrdered Regular Intervals crs: WellKnownText
+and 4 layers:
+ :tmin Float32 dims: X, Y (177×119)
+ :tmax Float32 dims: X, Y (177×119)
+ :prec Int16 dims: X, Y (177×119)
+ :wind Float32 dims: X, Y (177×119)
+
And plot scandinavia, with all borders included:
p = plot(scandinavia)
+borders!(p, denmark_border)
+borders!(p, norway_border)
+borders!(p, sweden_border)
+
+And save to netcdf - a single multi-layered file, and tif, which will write a file for each stack layer.
write("scandinavia.nc", scandinavia)
+write("scandinavia.tif", scandinavia)
(tmin = "scandinavia_tmin.tif", tmax = "scandinavia_tmax.tif", prec = "scandinavia_prec.tif", wind = "scandinavia_wind.tif")
Rasters.jl provides a range of other methods that are being added to over time. Where applicable these methods read and write lazily to and from disk-based arrays of common raster file types. These methods also work for entire RasterStacks
and RasterSeries
using the same syntax.
Plotting in Makie
2-D rasters in Makie
Plotting in Makie works somewhat differently than Plots, since the recipe system is different. You can pass a 2-D raster to any surface-like function (heatmap
, contour
, contourf
, or even surface
for a 3D plot) with ease.
using CairoMakie, Makie
+using Rasters, RasterDataSources, ArchGDAL
+A = Raster(WorldClim{BioClim}, 5) # this is a 3D raster, so is not accepted.
+B = A[:, :, 1] # this converts to a 2D raster which Makie accepts!
+figure = Figure()
+plot(figure[1, 1], B)
+contour(figure[1, 2], B)
+ax = Axis(figure[2, 1]; aspect = DataAspect())
+contourf!(ax, B)
+surface(figure[2, 2], B) # even a 3D plot works!
+figure
3-D rasters and RasterStacks in Makie
This interface is experimental, and unexported for that reason. It may break at any time!
Just as in Plots, 3D rasters are treated as a series of 2D rasters, which are tiled and plotted.
You can use Rasters.rplot
to visualize 3D rasters or RasterStacks in this way. An example is below:
stack = RasterStack(WorldClim{Climate}; month = 1)
+Rasters.rplot(stack; Axis = (aspect = DataAspect(),))
You can pass any theming keywords in, which are interpreted by Makie appropriately.
The plots seem a little squished here. We provide a Makie theme which makes text a little smaller and has some other space-efficient attributes:
CairoMakie.set_theme!(Rasters.theme_rasters())
+Rasters.rplot(stack)
Plotting with Observable
s
Rasters.rplot
should support Observable input out of the box, but the dimensions of that input must remain the same - i.e., the element names of a RasterStack must remain the same.
stack_obs = Observable(stack)
+fig = Rasters.rplot(stack_obs) # `stack` is the WorldClim climate data for January
+record(fig, "rplot.mp4", 1:12; framerate = 3) do i
+ stack_obs[] = RasterStack(WorldClim{Climate}; month = i)
+end
"rplot.mp4"
Rasters.rplot
— FunctionRasters.rplot([position::GridPosition], raster; kw...)
raster
may be a Raster
(of 2 or 3 dimensions) or a RasterStack
whose underlying rasters are 2 dimensional, or 3-dimensional with a singleton (length-1) third dimension.
Keywords
plottype = Makie.Heatmap
: The type of plot. Can be any Makie plot type which accepts aRaster
; in practice,Heatmap
,Contour
,Contourf
andSurface
are the best bets.axistype = Makie.Axis
: The type of axis. This can be anAxis
,Axis3
,LScene
, or even aGeoAxis
from GeoMakie.jl.X=X
: The X dimension of the raster.Y=Y
: The Y dimension of the raster.draw_colorbar = true
: Whether to draw a colorbar for the axis or not.colorbar_position = Makie.Right()
: Indicates which side of the axis the colorbar should be placed on. Can beMakie.Top()
,Makie.Bottom()
,Makie.Left()
, orMakie.Right()
.colorbar_padding = Makie.automatic
: The amound of padding between the colorbar and its axis. Ifautomatic
, then this is set to the width of the colorbar.title = Makie.automatic
: The titles of each plot. Ifautomatic
, these are set to the name of the band.xlabel = Makie.automatic
: The x-label for the axis. Ifautomatic
, set to the dimension name of the X-dimension of the raster.ylabel = Makie.automatic
: The y-label for the axis. Ifautomatic
, set to the dimension name of the Y-dimension of the raster.colorbarlabel = ""
: Usually nothing, but here if you need it. Sets the label on the colorbar.colormap = nothing
: The colormap for the heatmap. This can be set to a vector of colormaps (symbols, strings,cgrad
s) if plotting a 3D raster or RasterStack.colorrange = Makie.automatic
: The colormap for the heatmap. This can be set to a vector of(low, high)
if plotting a 3D raster or RasterStack.nan_color = :transparent
: The color whichNaN
values should take. Default to transparent.
Objects
Raster
Spatial raster data is essentially just an Array
. But Raster
wrappers allow treating them as an array that maintains its spatial index, crs and other metadata through all transformations. This means they can always be plotted and written to disk after applying most base Julia methods, and most broadcast
s.
Rasters.AbstractRaster
— TypeAbstractRaster <: DimensionalData.AbstractDimArray
Abstract supertype for objects that wrap an array (or location of an array) and metadata about its contents. It may be memory or hold a FileArray
, which holds the filename, and is only opened when required.
AbstractRaster
s inherit from AbstractDimArray
from DimensionalData.jl. They can be indexed as regular Julia arrays or with DimensionalData.jl Dimension
s. They will plot as a heatmap in Plots.jl with correct coordinates and labels, even after slicing with getindex
or view
. getindex
on a AbstractRaster
will always return a memory-backed Raster
.
Rasters.Raster
— TypeRaster <: AbsractRaster
+
+Raster(filepath::AbstractString, dims; kw...)
+Raster(A::AbstractArray{T,N}, dims; kw...)
+Raster(A::AbstractRaster; kw...)
A generic AbstractRaster
for spatial/raster array data. It may hold memory-backed arrays or FileArray
, that simply holds the String
path to an unopened file. This will only be opened lazily when it is indexed with getindex
or when read(A)
is called. Broadcasting, taking a view, reversing and most other methods do not load data from disk: they are applied later, lazily.
Keywords
dims
:Tuple
ofDimension
s for the array.lazy
: ABool
specifying if to load the stack lazily from disk.false
by default.name
:Symbol
name for the array, which will also retreive named layers ifRaster
is used on a multi-layered file like a NetCDF.missingval
: value reprsenting missing data, normally detected form the file. Set manually when you know the value is not specified or is incorrect. This will not change any values in the raster, it simply assigns which value is treated as missing. To replace all of the missing values in the raster, usereplace_missing
.metadata
:ArrayMetadata
object for the array, orNoMetadata()
.crs
: the coordinate reference system of the objectsXDim
/YDim
dimensions. Only set this if you know the detected crs is incrorrect, or it is not present in the file. Thecrs
is expected to be a GeoFormatTypes.jlCRS
orMixed
GeoFormat
type.mappedcrs
: the mapped coordinate reference system of the objectsXDim
/YDim
dimensions. forMapped
lookups these are the actual values of the index. ForProjected
lookups this can be used to index in eg.EPSG(4326)
lat/lon values, having it converted automatically. Only set this if the detectedmappedcrs
in incorrect, or the file does not have amappedcrs
, e.g. a tiff. Themappedcrs
is expected to be a GeoFormatTypes.jlCRS
orMixed
GeoFormat
type.dropband
: drop single band dimensions.true
by default.
Internal Keywords
In some cases it is possible to set these keywords as well.
data
: can replace the data in anAbstractRaster
refdims
:Tuple of
positionDimension
s the array was sliced from, defaulting to()
.
Rasters.Raster
— MethodRaster(T::Type{<:RasterDataSource}, [layer]; kw...) => Raster
Load a RasterDataSource
as an Raster
. T
and layers
are are passed to RasterDataSources.getraster
, while kw
args are for both getraster
and Raster
.
Keywords
month
: anInt
between1
and12
, usually forClimate
datasetsdate
: aDateTime
object, usually forWeather
datasets.res
: aString
resolution, for datasets with multiple resolutions.
MODIS
datasets require a specific set of keyword arguments:
lat
andlon
: Coordinates in decimal degrees (Float
s ) of the center of the rasterkm_ab
andkm_lr
: Kilometers above/below and left/right (Integer
s up to 100) of the center of the raster
Other Raster
keywords are passed to the Raster
constructor.
See the docs for RasterDatasources.getraster
for more specific details about data sources, layers and keyword arguments.
RasterStack
Spatial data often comes as a bundle of multiple named arrays, as in netcdf. RasterStack
can represent this, or multiple files organised in a similar way.
Rasters.AbstractRasterStack
— TypeAbstractRasterStack
Abstract supertype for objects that hold multiple AbstractRaster
s that share spatial dimensions.
They are NamedTuple
-like structures that may either contain NamedTuple
of AbstractRaster
s, string paths that will load AbstractRaster
s, or a single path that points to a file containing multiple layers, like NetCDF or HDF5. Use and syntax is similar or identical for all cases.
AbstractRasterStack
can hold layers that share some or all of their dimensions. They cannot have the same dimension with different length or spatial extent as another layer.
getindex
on an AbstractRasterStack
generally returns a memory backed standard Raster
. raster[:somelayer] |> plot
plots the layers array, while raster[:somelayer, X(1:100), Band(2)] |> plot
will plot the subset without loading the whole array.
getindex
on an AbstractRasterStack
with a key returns another stack with getindex
applied to all the arrays in the stack.
Rasters.RasterStack
— TypeRasterStack <: AbstrackRasterStack
+
+RasterStack(data...; name, kw...)
+RasterStack(data::Union{Vector,Tuple}; name, kw...)
+RasterStack(data::NamedTuple; kw...))
+RasterStack(s::AbstractRasterStack; kw...)
+RasterStack(s::AbstractRaster; layersfrom=Band, kw...)
+RasterStack(filename::AbstractString; kw...)
Load a file path or a NamedTuple
of paths as a RasterStack
, or convert arguments, a Vector
or NamedTuple
of Raster
s to RasterStack
.
Arguments
data
: ANamedTuple
ofRaster
s, or aVector
,Tuple
or splatted arguments ofRaster
. The latter options must pass aname
keyword argument.filename
: A file (such as netcdf or tif) to be loaded as a stack, or a directory path containing multiple files.
Keywords
name
: Used as stack layer names when aTuple
,Vector
or splat ofRaster
is passed in. Has no effect whenNameTuple
is used - theNamedTuple
keys are the layer names.metadata
: ADict
orDimensionalData.Metadata
object.refdims
:Tuple
ofDimension
that the stack was sliced from.layersfrom
:Dimension
to source stack layers from if the file is not already multi-layered.nothing
is default, so that a singleRasterStack(raster)
is a single layered stack.RasterStack(raster; layersfrom=Band)
will use the bands as layers.lazy
: ABool
specifying whether to load the stack lazily from disk.false
by default.dropband
: drop single band dimensions when creating stacks from filenames.true
by default.
files = (temp="temp.tif", pressure="pressure.tif", relhum="relhum.tif")
+stack = RasterStack(files; mappedcrs=EPSG(4326))
+stack[:relhum][Lat(Contains(-37), Lon(Contains(144))
Rasters.RasterStack
— MethodRasterStack(T::Type{<:RasterDataSource}, [layers::Union{Symbol,AbstractArray,Tuple}]; kw...) => RasterStack
Load a RasterDataSource
as an RasterStack
. T
and layers
are passed to RasterDataSources.getraster
, while kw
args are for both getraster
and RasterStack
.
Keywords
month
: anInt
between1
and12
, usually forClimate
datasets.date
: aDateTime
object, usually forWeather
datasets.res
: aString
resolution, for datasets with multiple resolutions.
MODIS
datasets require a specific set of keyword arguments:
lat
andlon
: Coordinates in decimal degrees (Float
s ) of the center of the rasterkm_ab
andkm_lr
: Kilometers above/below and left/right (Integer
s up to 100) of the center of the raster
Other RasterStack
keywords are passed to the RasterStack
constructor.
See the docs for RasterDatasources.getraster
for more specific details about data sources, layers and keyword arguments.
RasterSeries
A series is a meta-array that holds other files/data that is distributed over some dimension, often time. These files/data can be Raster
s or RasterStack
s.
Rasters.AbstractRasterSeries
— TypeAbstractRasterSeries <: DimensionalData.AbstractDimensionalArray
Abstract supertype for high-level DimensionalArray
that hold RasterStacks
, Rasters
, or the paths they can be loaded from. RasterSeries
are indexed with dimensions as with a AbstractRaster
. This is useful when you have multiple files containing rasters or stacks of rasters spread over dimensions like time and elevation.
As much as possible, implementations should facilitate loading entire directories and detecting the dimensions from metadata.
This allows syntax like below for a series of stacks of arrays:
RasterSeries[Time(Near(DateTime(2001, 1))][:temp][Y(Between(70, 150)), X(Between(-20,20))] |> plot`
RasterSeries
is the concrete implementation.
Rasters.RasterSeries
— TypeRasterSeries <: AbstractRasterSeries
+
+RasterSeries(rasters::AbstractArray{<:AbstractRaster}, dims; [refdims])
+RasterSeries(stacks::AbstractArray{<:AbstractRasterStack}, dims; [refdims])
+
+RasterSeries(paths::AbstractArray{<:AbstractString}, dims; child, duplicate_first, kw...)
+RasterSeries(path:::AbstractString, dims; ext, separator, child, duplicate_first, kw...)
Concrete implementation of AbstractRasterSeries
.
A RasterSeries
is an array of Raster
s or RasterStack
s, along some dimension(s).
Existing Raster
RasterStack
can be wrapped in a RasterSeries
, or new files can be loaded from an array of String
or from a single String
.
A single String
can refer to a whole directory, or the name of a series of files in a directory, sharing a common stem. The differnce between the filenames can be used as the lookup for the series.
For example, with some tifs at these paths :
"series_dir/myseries_2001-01-01T00:00:00.tif"
+"series_dir/myseries_2002-01-01T00:00:00.tif"
We can load a RasterSeries
with a DateTime
lookup:
julia> ser = RasterSeries("series_dir/myseries.tif", Ti(DateTime))
+2-element RasterSeries{Raster,1} with dimensions:
+ Ti Sampled{DateTime} DateTime[DateTime("2001-01-01T00:00:00"), DateTime("2002-01-01T00:00:00")] ForwardOrdered Irregular Points
The DateTime
suffix is parsed from the filenames. Using Ti(Int)
would try to parse integers intead.
Just using the directory will also work, unless there are other files mixed in it:
julia> ser = RasterSeries("series_dir", Ti(DateTime))
+2-element RasterSeries{Raster,1} with dimensions:
+ Ti Sampled{DateTime} DateTime[DateTime("2001-01-01T00:00:00"), DateTime("2002-01-01T00:00:00")] ForwardOrdered Irregular Points
Arguments
dims
: series dimension/s.
Keywords
When loading a series from a Vector of String
paths or a single String
path:
child
: constructor of child objects for use when filenames are passed in, can beRaster
orRasterStack
. Defaults toRaster
.duplicate_first::Bool
: wether to duplicate the dimensions and metadata of the first file with all other files. This can save load time with a large series where dimensions are identical.false
by default.lazy
: load files lazily,false
by default.kw
: keywords passed to the child constructorRaster
orRasterStack
.
When loading a series from a single String
path:
ext
: filename extension such as ".tiff" or ".nc". Use to specify a subset of files if only a directory path is passed in.separator
: separator used to split lookup elements from the rest of a filename. '_' by default.
Others:
refdims
: existing reference dimension/s, normally not required.
Rasters.RasterSeries
— MethodRasterSeries(T::Type{<:RasterDataSource}, [layers::Union{Symbol,AbstractArray,Tuple}]; kw...) => AbstractRasterSeries
Load a RasterDataSource
as an AbstractRasterSeries
. T
, args
are are passed to RasterDataSource.getraster
, while kw
args are for both getraster
and RasterSeries
.
Keywords
month
: aVector
or range ofInt
between1
and12
, usually forClimate
datasets.date
: aVector
ofDateTime
objects, usually forWeather
datasets.res
: aString
resolution, for datasets with multiple resolutions.
MODIS
datasets require a specific set of keyword arguments:
lat
andlon
: Coordinates in decimal degrees (Float
s ) of the center of the rasterkm_ab
andkm_lr
: Kilometers above/below and left/right (Integer
s up to 100) of the center of the raster
Other RasterSeries
keywords are passed to the RasterSeries
constructor.
See the docs for RasterDatasources.getraster
for more specific details about data sources, layers and keyword arguments.
Dimensions
Rasters uses X
, Y
, and Z
dimensions from DimensionalData.jl to represent spatial directions like longitude, latitude and the vertical dimension, and subset data with them. Ti
is used for time, and Band
represent bands. Other dimensions can have arbitrary names, but will be treated generically. See DimensionalData.jl for more details on how they work.
Rasters.Band
— TypeBand <: Dimension
+
+Band(val=:)
Band Dimension
for multi-band rasters.
Example:
banddim = Band(10:10:100)
+# Or
+val = A[Band(1)]
+# Or
+mean(A; dims=Band)
Lookup Arrays
These specify properties of the index associated with e.g. the X and Y dimension. Rasters.jl defines additional lookup arrays: Projected
to handle dimensions with projections, and Mapped
where the projection is mapped to another projection like EPSG(4326)
. Mapped
is largely designed to handle NetCDF dimensions, especially with Explicit
spans.
Rasters.AbstractProjected
— TypeAbstractProjected <: AbstractSampled
Abstract supertype for projected index lookups.
Rasters.Projected
— TypeProjected <: AbstractProjected
+
+Projected(order, span, sampling, crs, mappedcrs)
+Projected(; order=AutoOrder(), span=AutoSpan(), sampling=AutoSampling(), crs, mappedcrs=nothing)
An AbstractSampled
LookupArray
with projections attached.
Fields and behaviours are identical to Sampled
with the addition of crs
and mappedcrs
fields.
If both crs
and mappedcrs
fields contain CRS data (in a GeoFormat
wrapper from GeoFormatTypes.jl) the selector inputs and plot axes will be converted from and to the specified mappedcrs
projection automatically. A common use case would be to pass mappedcrs=EPSG(4326)
to the constructor when loading eg. a GDALarray:
GDALarray(filename; mappedcrs=EPSG(4326))
The underlying crs
will be detected by GDAL.
If mappedcrs
is not supplied (ie. mappedcrs=nothing
), the base index will be shown on plots, and selectors will need to use whatever format it is in.
Rasters.Mapped
— TypeMapped <: AbstractProjected
+
+Mapped(order, span, sampling, crs, mappedcrs)
+Mapped(; order=AutoOrder(), span=AutoSpan(), sampling=AutoSampling(), crs=nothing, mappedcrs)
An AbstractSampled
LookupArray
, where the dimension index has been mapped to another projection, usually lat/lon or EPSG(4326)
. Mapped
matches the dimension format commonly used in netcdf files.
Fields and behaviours are identical to Sampled
with the addition of crs
and mappedcrs
fields.
The mapped dimension index will be used as for Sampled
, but to save in another format the underlying crs
may be used to convert it.
Data sources
Rasters.jl uses a number of backends to load raster data. Raster
, RasterStack
and RasterSeries
will detect which backend to use for you, automatically.
GRD
R GRD files can be loaded natively, using Julias MMap
- which means they are very fast, but are not compressed. They are always 3 dimensional, and have Y
, X
and Band
dimensions.
NetCDF
NetCDF .nc
files are loaded using NCDatasets.jl. Layers from files can be loaded as Raster("filename.nc"; key=:layername)
. Without key
the first layer is used. RasterStack("filename.nc")
will use all netcdf variables in the file that are not dimensions as layers.
NetCDF layers can have arbitrary dimensions. Known, common dimension names are converted to X
, Y
Z
, and Ti
, otherwise Dim{:layername}
is used. Layers in the same file may also have different dimensions.
NetCDF files still have issues loading directly from disk for some operations. Using read(ncstack)
may help.
GDAL
All files GDAL can access, such as .tiff
and .asc
files, can be loaded, using ArchGDAL.jl. These are generally best loaded as Raster("filename.tif")
, but can be loaded as RasterStack("filename.tif"; layersfrom=Band)
, taking layers from the Band
dimension, which is also the default.
SMAP
The Soil Moisture Active-Passive dataset provides global layers of soil moisture, temperature and other related data, in a custom HDF5 format. Layers are always 2 dimensional, with Y
and X
dimensions.
These can be loaded as multi-layered RasterStack("filename.h5")
. Individual layers can be loaded as Raster("filename.h5"; key=:layerkey)
, without key
the first layer is used.
Rasters.smapseries
— Functionsmapseries(filenames::AbstractString; kw...)
+smapseries(filenames::Vector{<:AbstractString}, dims=nothing; kw...)
RasterSeries
loader for SMAP files and whole folders of files, organised along the time dimension. Returns a RasterSeries
.
Arguments
filenames
: AString
path to a directory of SMAP files, or a vector ofString
paths to specific files.dims
:Tuple
containingTi
dimension for the series. Automatically generated formfilenames
unless passed in.
Keywords
kw
: Passed toRasterSeries
.
Writing file formats to disk
Files can be written to disk in all formats other than SMAP HDF5 using write("filename.ext", A)
. See the docs for write
. They can (with some caveats) be written to different formats than they were loaded in as, providing file-type conversion for spatial data.
Some metadata may be lost in formats that store little metadata, or where metadata conversion has not been completely implemented.
RasterDataSources.jl integration
RasterDataSources.jl standardises the download of common raster data sources, with a focus on datasets used in ecology and the environmental sciences. RasterDataSources.jl is tightly integrated into Rasters.jl, so that datsets and keywords can be used directly to download and load data as a Raster
, RasterStack
, or RasterSeries
.
using Rasters, RasterDataSources, ArchGDAL, Plots, Dates
+A = Raster(WorldClim{Climate}, :tavg; month=June)
+plot(A)
+
+See the docs for Raster
, RasterStack
and RasterSeries
, and the docs for RasterDataSources.getraster
for syntax to specify various data sources.
Exported functions
Rasters.jl is a direct extension of DimensionalData.jl. See DimensionalData.jl docs for the majority of types and functions that can be used in Rasters.jl.
Functions more specific to geospatial data are included in Rasters.jl, and listed below.
Rasters.aggregate
— Functionaggregate(method, object, scale; filename, progress, skipmissing)
Aggregate a Raster
, or all arrays in a RasterStack
or RasterSeries
, by scale
using method
.
Arguments
method
: a function such asmean
orsum
that can combine the value of multiple cells to generate the aggregated cell, or aLocus
likeStart()
orCenter()
that specifies where to sample from in the interval.object
: Object to aggregate, likeAbstractRasterSeries
,AbstractStack
,AbstractRaster
orDimension
.scale
: the aggregation factor, which can be an integer, a tuple of integers for each dimension, or anyDimension
,Selector
orInt
combination you can usually use ingetindex
. Using aSelector
will determine the scale by the distance from the start of the index.
When the aggregation scale
of is larger than the array axis, the length of the axis is used.
Keywords
skipmissingval
: iftrue
, anymissingval
will be skipped during aggregation, so that only areas of all missing values will be aggregated tomissingval
. Iffalse
, any aggegrated area containing amissingval
will be assignedmissingval
.filename
: a filename to write to directly, useful for large files.suffix
: a string or value to append to the filename. A tuple ofsuffix
will be applied to stack layers.keys(stack)
are the default.progress
: show a progress bar,true
by default,false
to hide.
Example
using Rasters, RasterDataSources, Statistics, Plots
+using Rasters: Center
+st = read(RasterStack(WorldClim{Climate}; month=1))
+ag = aggregate(Center(), st, (Y(20), X(20)); skipmissingval=true, progress=false)
+plot(ag)
+savefig("build/aggregate_example.png"); nothing
+# output
+
Note: currently it is faster to aggregate over memory-backed arrays. Use read
on src
before use where required.
Rasters.aggregate!
— Functionaggregate!(method, dst::AbstractRaster, src::AbstractRaster, scale; skipmissingval=false)
Aggregate array src
to array dst
by scale
, using method
.
Arguments
method
: a function such asmean
orsum
that can combine the value of multiple cells to generate the aggregated cell, or aLocus
likeStart()
orCenter()
that species where to sample from in the interval.scale
: the aggregation factor, which can be an integer, a tuple of integers for each dimension, or anyDimension
,Selector
orInt
combination you can usually use ingetindex
. Using aSelector
will determine the scale by the distance from the start of the index in thesrc
array.
When the aggregation scale
of is larger than the array axis, the length of the axis is used.
Keywords
progress
: show a progress bar.skipmissingval
: iftrue
, anymissingval
will be skipped during aggregation, so that only areas of all missing values will be aggregated tomissingval
. Iffalse
, any aggegrated area containing amissingval
will be assignedmissingval
.
Note: currently it is much faster to aggregate over memory-backed arrays. Use read
on src
before use where required.
Rasters.boolmask
— Functionboolmask(obj::Raster; [missingval])
+boolmask(obj; [to, res, size])
Create a mask array of Bool
values, from another Raster
. An AbstractRasterStack
or AbstractRasterSeries
are also accepted, but a mask is taken of the first layer or object not all of them.
The array returned from calling boolmask
on a AbstractRaster
is a Raster
with the same dimensions as the original array and a missingval
of false
.
Arguments
obj
: aRaster
, a GeoInterface.jl geometry, or a vector or table of geometries.
Raster
/ RasterStack
Keywords
missingval
: The missing value of the source array, with defaultmissingval(raster)
.
Geometry keywords
to
: aRaster
,RasterStack
,Tuple
ofDimension
orExtents.Extent
. If noto
object is provided the extent will be calculated from the geometries, Additionally, when noto
object or anExtent
is passed forto
, thesize
orres
keyword must also be used.res
: the resolution of the dimensions, aReal
orTuple{<:Real,<:Real}
. Only required whento
is not used or is anExtents.Extent
, andsize
is not used.size
: the size of the output array, as aTuple{Int,Int}
or singleInt
for a square. Only required whento
is not used or is anExtents.Extent
, andres
is not used.crs
: acrs
which will be attached to the resulting raster whento
not passed or is anExtent
. Otherwise the crs fromto
is used.shape
: Forcedata
to be treated as:polygon
,:line
or:point
geometries. using points or lines as polygons may have unexpected results.boundary
: for polygons, include pixels where the:center
is inside the polygon, where the polygon:touches
the pixel, or that are completely:inside
the polygon. The default is:center
.
And specifically for shape=:polygon
:
boundary
: include pixels where the:center
is inside the polygon, where the line:touches
the pixel, or that are completely:inside
inside the polygon. The default is:center
.
For tabular data, feature collections and other iterables
collapse
: iftrue
, collapse all geometry masks into a single mask. Otherwise return a Raster with an additionalgeometry
dimension, so that each slice along this axis is the mask of thegeometry
opbject of each row of the table, feature in the feature collection, or just each geometry in the iterable.
Example
using Rasters, RasterDataSources, ArchGDAL, Plots, Dates
+wc = Raster(WorldClim{Climate}, :prec; month=1)
+boolmask(wc) |> plot
+
+savefig("build/boolmask_example.png"); nothing
+
+# output
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.classify
— Functionclassify(x, pairs; lower=(>=), upper=(<), others=nothing)
+classify(x, pairs...; lower, upper, others)
Create a new array with values in x
classified by the values in pairs
.
pairs
can hold tuples fo values (2, 3)
, a Fix2
function e.g. <=(1)
, a Tuple
of Fix2
e.g. (>=(4), <(7))
, or an IntervalSets.jl interval, e.g. 3..9
or OpenInterval(10, 12)
. pairs
can also be a n * 3
matrix where each row is lower bounds, upper bounds, replacement.
If tuples or a Matrix
are used, the lower
and upper
keywords define how the lower and upper boundaries are chosen.
If others
is set other values not covered in pairs
will be set to that values.
Arguments
x
: aRaster
orRasterStack
pairs
: each pair contains a value and a replacement, a tuple of lower and upper range and a replacement, or a Tuple ofFix2
like(>(x), <(y)
.
Keywords
lower
: Which comparison (<
or<=
) to use for lower values, ifFix2
are not used.upper
: Which comparison (>
or>=
) to use for upper values, ifFix2
are not used.others
: A value to assign to all values not included inpairs
. Passingnothing
(the default) will leave them unchanged.
Example
using Rasters, RasterDataSources, ArchGDAL, Plots
+A = Raster(WorldClim{Climate}, :tavg; month=1)
+classes = <=(15) => 10,
+ 15..25 => 20,
+ 25..35 => 30,
+ >(35) => 40
+classified = classify(A, classes; others=0, missingval=0)
+plot(classified; c=:magma)
+
+savefig("build/classify_example.png"); nothing
+
+# output
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.classify!
— Functionclassify!(x, pairs...; lower, upper, others)
+classify!(x, pairs; lower, upper, others)
Classify the values of x
in-place, by the values in pairs
.
If Fix2
is not used, the lower
and upper
keywords
If others
is set other values not covered in pairs
will be set to that values.
Arguments
x
: aRaster
orRasterStack
pairs
: each pair contains a value and a replacement, a tuple of lower and upper range and a replacement, or a Tuple ofFix2
like(>(x), <(y)
.
Keywords
lower
: Which comparison (<
or<=
) to use for lower values, ifFix2
are not used.upper
: Which comparison (>
or>=
) to use for upper values, ifFix2
are not used.others
: A value to assign to all values not included inpairs
. Passingnothing
(the default) will leave them unchanged.
Example
classify!
to disk, with key steps:
- copying a tempory file so we don't write over the RasterDataSources.jl version.
- use
open
withwrite=true
to open the file with disk-write permissions. - use
Float32
like10.0f0
for all our replacement values andother
, because the file is stored asFloat32
. Attempting to write some other type will fail.
using Rasters, RasterDataSources, ArchGDAL, Plots
+# Download and copy the file
+filename = getraster(WorldClim{Climate}, :tavg; month=6)
+tempfile = tempname() * ".tif"
+cp(filename, tempfile)
+# Define classes
+classes = (5, 15) => 10,
+ (15, 25) => 20,
+ (25, 35) => 30,
+ >=(35) => 40
+# Open the file with write permission
+open(Raster(tempfile); write=true) do A
+ classify!(A, classes; others=0)
+end
+# Open it again to plot the changes
+plot(Raster(tempfile); c=:magma)
+
+savefig("build/classify_bang_example.png"); nothing
+
+# output
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.coverage
— Functioncoverage(mode, geom; [to, res, size, scale, verbose, progress])
+coverage(geom; [to, mode, res, size, scale, verbose, progress])
Calculate the area of a raster covered by GeoInterface.jl compatible geomtry geom
, as a fraction.
Each pixel is assigned a grid of points (by default 10 x 10) that are each checked to be inside the geometry. The sum divided by the number of points to give coverage.
In pracice, most pixel coverage is not calculated this way - shortcuts that produce the same result are taken wherever possible.
If geom
is an AbstractVector
or table, the mode
keyword will determine how coverage is combined.
Keywords
mode
: method for combining multiple geometries -union
orsum
.union
(the default) gives the areas covered by all geometries. Usefull in spatial coverage where overlapping regions should not be counted twice. The returned raster will containFloat64
values between0.0
and1.0
.sum
gives the summed total of the areas covered by all geometries, as in taking the sum of runningcoverage
separately on all geometries. The returned values are positiveFloat64
.
For a single geometry, the
mode
keyword has no effect - the result is the same.scale
:Integer
scale of pixel subdivision. The default of10
means each pixel has 10 x 10 or 100 points that contribute to coverage. Using100
means 10,000 points contribute. Performance will decline asscale
increases. Memory use will grow byscale^2
whenmode=:union
.progress
: show a progress bar,true
by default,false
to hide.vebose
: whether to print messages about potential problems.true
by default.
to
: aRaster
,RasterStack
,Tuple
ofDimension
orExtents.Extent
. If noto
object is provided the extent will be calculated from the geometries, Additionally, when noto
object or anExtent
is passed forto
, thesize
orres
keyword must also be used.size
: the size of the output array, as aTuple{Int,Int}
or singleInt
for a square. Only required whento
is not used or is anExtents.Extent
, andres
is not used.res
: the resolution of the dimensions, aReal
orTuple{<:Real,<:Real}
. Only required whento
is not used or is anExtents.Extent
, andsize
is not used.
Rasters.coverage!
— Functioncoverage!(A, geom; [mode, scale])
Calculate the area of a raster covered by GeoInterface.jl compatible geomtry geom
, as a fraction.
Each pixel is assigned a grid of points (by default 10 x 10) that are each checked to be inside the geometry. The sum divided by the number of points to give coverage.
In pracice, most pixel coverage is not calculated this way - shortcuts that produce the same result are taken wherever possible.
If geom
is an AbstractVector
or table, the mode
keyword will determine how coverage is combined.
Keywords
mode
: method for combining multiple geometries -union
orsum
.union
(the default) gives the areas covered by all geometries. Usefull in spatial coverage where overlapping regions should not be counted twice. The returned raster will containFloat64
values between0.0
and1.0
.sum
gives the summed total of the areas covered by all geometries, as in taking the sum of runningcoverage
separately on all geometries. The returned values are positiveFloat64
.
For a single geometry, the
mode
keyword has no effect - the result is the same.scale
:Integer
scale of pixel subdivision. The default of10
means each pixel has 10 x 10 or 100 points that contribute to coverage. Using100
means 10,000 points contribute. Performance will decline asscale
increases. Memory use will grow byscale^2
whenmode=:union
.progress
: show a progress bar,true
by default,false
to hide.vebose
: whether to print messages about potential problems.true
by default.
Rasters.convertlookup
— Functionconvertlookup(dstlookup::Type{<:LookupArray}, x)
Convert the dimension lookup between Projected
and Mapped
. Other dimension lookups pass through unchanged.
This is used to e.g. save a netcdf file to GeoTiff.
Rasters.crop
— Functioncrop(x; to)
+crop(xs...; to)
Crop one or multiple AbstractRaster
or AbstractRasterStack
x
to match the size of the object to
, or smallest of any dimensions that are shared.
crop
is lazy, using a view
into the object rather than alocating new memory.
Keywords
to
: the object to crop to. Ifto
keyword is passed, the smallest shared area of allxs
is used.touches
:true
orfalse
. Whether to useTouches
wraper on the object extent. When lines need to be included in e.g. zonal statistics,true
shoudle be used.
As crop
is lazy, filename
and suffix
keywords are not used.
Example
Crop to another raster:
using Rasters, RasterDataSources, Plots
+evenness = Raster(EarthEnv{HabitatHeterogeneity}, :evenness)
+rnge = Raster(EarthEnv{HabitatHeterogeneity}, :range)
+
+# Roughly cut out New Zealand from the evenness raster
+nz_bounds = X(165 .. 180), Y(-50 .. -32)
+nz_evenness = evenness[nz_bounds...]
+
+# Crop range to match evenness
+nz_range = crop(rnge; to=nz_evenness)
+plot(nz_range)
+
+savefig("build/nz_crop_example.png")
+nothing
+
+# output
Crop to a polygon:
using Rasters, RasterDataSources, Plots, Dates, Shapefile, Downloads
+
+# Download a borders shapefile
+shapefile_url = "https://github.com/nvkelso/natural-earth-vector/raw/master/10m_cultural/ne_10m_admin_0_countries.shp"
+shapefile_name = "boundary.shp"
+isfile(shapefile_name) || Downloads.download(shapefile_url, shapefile_name)
+shp = Shapefile.Handle(shapefile_name).shapes[6]
+
+evenness = Raster(EarthEnv{HabitatHeterogeneity}, :evenness)
+argentina_evenness = crop(evenness; to=shp)
+plot(argentina_evenness)
+
+savefig("build/argentina_crop_example.png"); nothing
+
+# output
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.crs
— Functioncrs(x)
Get the projected coordinate reference system of a Y
or X
Dimension
, or of the Y
/X
dims of an AbstractRaster
.
For Mapped
lookup this may be nothing
as there may be no projected coordinate reference system at all. See setcrs
to set it manually.
Rasters.disaggregate
— Functiondisaggregate(method, object, scale; filename, progress, keys)
Disaggregate array, or all arrays in a stack or series, by some scale.
Arguments
method
: a function such asmean
orsum
that can combine the value of multiple cells to generate the aggregated cell, or aLocus
likeStart()
orCenter()
that species where to sample from in the interval.object
: Object to aggregate, likeAbstractRasterSeries
,AbstractStack
,AbstractRaster
or aDimension
.scale
: the aggregation factor, which can be an integer, a tuple of integers for each dimension, or anyDimension
,Selector
orInt
combination you can usually use ingetindex
. Using aSelector
will determine the scale by the distance from the start of the index.
Keywords
progress
: show a progress bar.
Note: currently it is faster to aggregate over memory-backed arrays. Use read
on src
before use where required.
Rasters.disaggregate!
— Functiondisaggregate!(method, dst::AbstractRaster, src::AbstractRaster, filename, scale)
Disaggregate array src
to array dst
by some scale, using method
.
method
: a function such asmean
orsum
that can combine the value of multiple cells to generate the aggregated cell, or aLocus
likeStart()
orCenter()
that species where to sample from in the interval.scale
: the aggregation factor, which can be an integer, a tuple of integers for each dimension, or anyDimension
,Selector
orInt
combination you can usually use ingetindex
. Using aSelector
will determine the scale by the distance from the start of the index in thesrc
array.
Note: currently it is faster to aggregate over memory-backed arrays. Use read
on src
before use where required.
Rasters.extend
— Functionextend(xs...; [to])
+extend(xs; [to])
+extend(x::Union{AbstractRaster,AbstractRasterStack}; to, kw...)
Extend one or multiple AbstractRaster
to match the area covered by all xs
, or by the keyword argument to
.
Keywords
to
: the Raster or dims to extend to. If noto
keyword is passed, the largest shared area of allxs
is used.touches
:true
orfalse
. Whether to useTouches
wraper on the object extent. When lines need to be included in e.g. zonal statistics,true
shoudle be used.filename
: a filename to write to directly, useful for large files.suffix
: a string or value to append to the filename. A tuple ofsuffix
will be applied to stack layers.keys(stack)
are the default.
using Rasters, RasterDataSources, Plots
+evenness = Raster(EarthEnv{HabitatHeterogeneity}, :evenness)
+rnge = Raster(EarthEnv{HabitatHeterogeneity}, :range)
+
+# Roughly cut out South America
+sa_bounds = X(-88 .. -32), Y(-57 .. 13)
+sa_evenness = evenness[sa_bounds...]
+
+# Extend range to match the whole-world raster
+sa_range = extend(sa_evenness; to=rnge)
+plot(sa_range)
+
+savefig("build/extend_example.png")
+nothing
+# output
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.extract
— Functionextract(x, geoms; atol)
Extracts the value of Raster
or RasterStack
at given points, returning an iterable of NamedTuple
with properties for :geometry
and raster or stack layer values.
Note that if objects have more dimensions than the length of the point tuples, sliced arrays or stacks will be returned instead of single values.
Arguments
x
: aRaster
orRasterStack
to extract values from.geoms
: GeoInterface.jl compatible geometries, or tables or iterables of geometries.
Keywords
atol
: a tolorerance for floating point lookup values for when theLookupArray
containsPoints
.atol
is ignored forIntervals
.
Example
Here we extact points matching the occurrence of the Mountain Pygmy Possum, Burramis parvus. This could be used to fit a species distribution model.
using Rasters, RasterDataSources, ArchGDAL, GBIF2, CSV
+
+# Get a stack of BioClim layers, and replace missing values with `missing`
+st = RasterStack(WorldClim{BioClim}, (1, 3, 5, 7, 12)) |> replace_missing
+
+# Download some occurrence data
+obs = GBIF2.occurrence_search("Burramys parvus"; limit=5, year="2009")
+
+# Convert observations to points
+pnts = collect((o.decimalLongitude, o.decimalLatitude) for o in obs if !ismissing(o.decimalLongitude))
+
+# use `extract` to get values for all layers at each observation point.
+# We `collect` to get a `Vector` from the lazy iterator.
+collect(extract(st, pnts))
+
+# output
+5-element Vector{NamedTuple{(:geometry, :bio1, :bio3, :bio5, :bio7, :bio12)}}:
+ (geometry = (0.21, 40.07), bio1 = 17.077084f0, bio3 = 41.20417f0, bio5 = 30.1f0, bio7 = 24.775f0, bio12 = 446.0f0)
+ (geometry = (0.03, 39.97), bio1 = 17.076923f0, bio3 = 39.7983f0, bio5 = 29.638462f0, bio7 = 24.153847f0, bio12 = 441.0f0)
+ (geometry = (0.03, 39.97), bio1 = 17.076923f0, bio3 = 39.7983f0, bio5 = 29.638462f0, bio7 = 24.153847f0, bio12 = 441.0f0)
+ (geometry = (0.52, 40.37), bio1 = missing, bio3 = missing, bio5 = missing, bio7 = missing, bio12 = missing)
+ (geometry = (0.32, 40.24), bio1 = 16.321388f0, bio3 = 41.659454f0, bio5 = 30.029825f0, bio7 = 25.544561f0, bio12 = 480.0f0)
Rasters.mappedcrs
— Functionmappedcrs(x)
Get the mapped coordinate reference system for the Y
/X
dims of an array.
In Projected
lookup this is used to convert Selector
values form the mappedcrs defined projection to the underlying projection, and to show plot axes in the mapped projection.
In Mapped
lookup this is the coordinate reference system of the index values. See setmappedcrs
to set it manually.
Rasters.mappedbounds
— Functionmappedbounds(x)
Get the bounds converted to the mappedcrs
value.
Whithout ArchGDAL loaded, this is just the regular bounds.
Rasters.mappedindex
— Functionmappedindex(x)
Get the index value of a dimension converted to the mappedcrs
value.
Whithout ArchGDAL loaded, this is just the regular dim value.
Rasters.mask
— Functionmask(A:AbstractRaster; with, missingval=missingval(A))
+mask(x; with)
Return a new array with values of A
masked by the missing values of with
, or by the shape of with
, if with
is a geometric object.
Arguments
x
: aRaster
orRasterStack
Keywords
with
: anAbstractRaster
, or any GeoInterface.jl compatible objects or table. The coordinate reference system of the point must matchcrs(A)
.missingval
: the missing value to use in the returned file.filename
: a filename to write to directly, useful for large files.suffix
: a string or value to append to the filename. A tuple ofsuffix
will be applied to stack layers.keys(stack)
are the default.
Geometry keywords
These can be used when with
is a GeoInterface.jl compatible object:
shape
: Forcedata
to be treated as:polygon
,:line
or:point
geometries. using points or lines as polygons may have unexpected results.boundary
: for polygons, include pixels where the:center
is inside the polygon, where the polygon:touches
the pixel, or that are completely:inside
the polygon. The default is:center
.
Example
Mask an unmasked AWAP layer with a masked WorldClim layer, by first resampling the mask.
using Rasters, RasterDataSources, ArchGDAL, Plots, Dates
+
+# Load and plot the file
+awap = read(Raster(AWAP, :tmax; date=DateTime(2001, 1, 1)))
+a = plot(awap; clims=(10, 45))
+
+# Create a mask my resampling a worldclim file
+wc = Raster(WorldClim{Climate}, :prec; month=1)
+wc_mask = resample(wc; to=awap)
+
+# Mask
+awap_masked = mask(awap; with=wc_mask)
+b = plot(awap_masked; clims=(10, 45))
+
+savefig(a, "build/mask_example_before.png");
+savefig(b, "build/mask_example_after.png"); nothing
+# output
+
Before mask
:
After mask
:
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.mask!
— Functionmask!(x; with, missingval=missingval(A))
Mask A
by the missing values of with
, or by all values outside with
if it is a polygon.
If with
is a polygon, creates a new array where points falling outside the polygon have been replaced by missingval(A)
.
Return a new array with values of A
masked by the missing values of with
, or by a polygon.
Arguments
x
: aRaster
orRasterStack
.
Keywords
with
: anotherAbstractRaster
, aAbstractVector
ofTuple
points, or any GeoInterface.jlAbstractGeometry
. The coordinate reference system of the point must matchcrs(A)
.missingval
: the missing value to write to A in masked areas, by defaultmissingval(A)
.
Example
Mask an unmasked AWAP layer with a masked WorldClim layer, by first resampling the mask to match the size and projection.
using Rasters, RasterDataSources, ArchGDAL, Plots, Dates
+
+# Load and plot the file
+awap = read(RasterStack(AWAP, (:tmin, :tmax); date=DateTime(2001, 1, 1)))
+a = plot(awap; clims=(10, 45), c=:imola)
+
+# Create a mask my resampling a worldclim file
+wc = Raster(WorldClim{Climate}, :prec; month=1)
+wc_mask = resample(wc; to=awap)
+
+# Mask
+mask!(awap; with=wc_mask)
+b = plot(awap; clims=(10, 45))
+
+savefig(a, "build/mask_bang_example_before.png");
+savefig(b, "build/mask_bang_example_after.png"); nothing
+
+# output
+
Before mask!
:
After mask!
:
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.missingval
— Functionmissingval(x)
Returns the value representing missing data in the dataset
Rasters.missingmask
— Functionmissingmask(obj::Raster; kw...)
+missingmask(obj; [to, res, size, collapse])
Create a mask array of missing
and true
values, from another Raster
. AbstractRasterStack
or AbstractRasterSeries
are also accepted, but a mask is taken of the first layer or object not all of them.
For AbstractRaster
the default missingval
is missingval(A)
, but others can be chosen manually.
The array returned from calling missingmask
on a AbstractRaster
is a Raster
with the same size and fields as the original array.
Keywords
to
: aRaster
,RasterStack
,Tuple
ofDimension
orExtents.Extent
. If noto
object is provided the extent will be calculated from the geometries, Additionally, when noto
object or anExtent
is passed forto
, thesize
orres
keyword must also be used.res
: the resolution of the dimensions, aReal
orTuple{<:Real,<:Real}
. Only required whento
is not used or is anExtents.Extent
, andsize
is not used.size
: the size of the output array, as aTuple{Int,Int}
or singleInt
for a square. Only required whento
is not used or is anExtents.Extent
, andres
is not used.crs
: acrs
which will be attached to the resulting raster whento
not passed or is anExtent
. Otherwise the crs fromto
is used.shape
: Forcedata
to be treated as:polygon
,:line
or:point
geometries. using points or lines as polygons may have unexpected results.boundary
: for polygons, include pixels where the:center
is inside the polygon, where the polygon:touches
the pixel, or that are completely:inside
the polygon. The default is:center
.
Example
using Rasters, RasterDataSources, ArchGDAL, Plots, Dates
+wc = Raster(WorldClim{Climate}, :prec; month=1)
+missingmask(wc) |> plot
+
+savefig("build/missingmask_example.png"); nothing
+
+# output
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.mosaic
— Functionmosaic(f, regions...; missingval, atol)
+mosaic(f, regions; missingval, atol)
Combine regions
into a single raster.
Arguments
f
: A reducing function (mean
,sum
,first
,last
etc.) for values whereregions
overlap.regions
: Iterable or splattedRaster
orRasterStack
.
Keywords
missingval
: Fills empty areas, and defualts to themissingval
of the first region.atol
: Absolute tolerance for comparison between index values. This is often required due to minor differences in range values due to floating point error. It is not applied to non-float dimensions. A tuple of tolerances may be passed, matching the dimension order.filename
: a filename to write to directly, useful for large files.suffix
: a string or value to append to the filename. A tuple ofsuffix
will be applied to stack layers.keys(stack)
are the default.
If your mosaic has has apparent line errors, increase the atol
value.
Example
Here we cut out Australia and Africa from a stack, and join them with mosaic
.
using Rasters, RasterDataSources, ArchGDAL, Plots
+st = RasterStack(WorldClim{Climate}; month=1);
+
+africa = st[X(-20.0 .. 60.0), Y(-40.0 .. 35.0)]
+a = plot(africa)
+
+aus = st[X(100.0 .. 160.0), Y(-50.0 .. -10.0)]
+b = plot(aus)
+
+# Combine with mosaic
+mos = mosaic(first, aus, africa)
+c = plot(mos)
+
+savefig(a, "build/mosaic_example_africa.png")
+savefig(b, "build/mosaic_example_aus.png")
+savefig(c, "build/mosaic_example_combined.png")
+nothing
+# output
+
Individual continents
Mosaic of continents
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.mosaic!
— Functionmosaic!(f, x, regions...; missingval, atol)
+mosaic!(f, x, regions::Tuple; missingval, atol)
Combine regions
in x
using the function f
.
Arguments
f
a function (e.g.mean
,sum
,first
orlast
) that is applied to values whereregions
overlap.x
: ARaster
orRasterStack
. May be a an opened disk-basedRaster
, the result will be written to disk. With the current algorithm, the read speed is slow.regions
: source objects to be joined. These should be memory-backed (useread
first), or may experience poor performance. If all objects have the same extent,mosaic
is simply a merge.
Keywords
missingval
: Fills empty areas, and defualts to the `missingval/ of the first layer.atol
: Absolute tolerance for comparison between index values. This is often required due to minor differences in range values due to floating point error. It is not applied to non-float dimensions. A tuple of tolerances may be passed, matching the dimension order.
Example
Cut out Australia and Africa stacks, then combined them into a single stack.
using Rasters, RasterDataSources, ArchGDAL, Statistics, Plots
+st = read(RasterStack(WorldClim{Climate}; month=1))
+aus = st[X=100.0 .. 160.0, Y=-50.0 .. -10.0]
+africa = st[X=-20.0 .. 60.0, Y=-40.0 .. 35.0]
+mosaic!(first, st, aus, africa)
+plot(st)
+savefig("build/mosaic_bang_example.png")
+nothing
+# output
+
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.points
— Functionpoints(A::AbstractRaster; dims=(YDim, XDim), ignore_missing) => Array{Tuple}
Returns a generator of the points in A
for dimensions in dims
, where points are a tuple of the values in each specified dimension index.
Keywords
dims
the dimensions to return points from. The first slice of other layers will be used.ignore_missing
: wether to ignore missing values in the array when considering points. Iftrue
, all points in the dimensions will be returned, iffalse
only the points that are not=== missingval(A)
will be returned.
The order of dims
determines the order of the points.
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.rasterize
— Functionrasterize([reducer], data; kw...)
Rasterize a GeoInterface.jl compatable geometry or feature, or a Tables.jl table with a :geometry
column of GeoInterface.jl objects, or X
, Y
points columns.
Arguments
reducer
: a reducing function to reduce the fill value for all geometries that cover or touch a pixel down to a single value. The default islast
. Any that takes an iterable and returns a single value will work, including custom functions. However, there are optimisations for built-in methods includingsum
,first
,last
,minimum
,maximum
,extrema
andStatistics.mean
. These may be an order of magnitude or more faster thancount
is a special-cased as it does not need a fill value.data
: a GeoInterface.jlAbstractGeometry
, or a nestedVector
ofAbstractGeometry
, or a Tables.jl compatible object containing a:geometry
column or points and values columns.
Keywords
These are detected automatically from data
where possible.
to
: aRaster
,RasterStack
,Tuple
ofDimension
orExtents.Extent
. If noto
object is provided the extent will be calculated from the geometries, Additionally, when noto
object or anExtent
is passed forto
, thesize
orres
keyword must also be used.res
: the resolution of the dimensions, aReal
orTuple{<:Real,<:Real}
. Only required whento
is not used or is anExtents.Extent
, andsize
is not used.size
: the size of the output array, as aTuple{Int,Int}
or singleInt
for a square. Only required whento
is not used or is anExtents.Extent
, andres
is not used.crs
: acrs
which will be attached to the resulting raster whento
not passed or is anExtent
. Otherwise the crs fromto
is used.shape
: Forcedata
to be treated as:polygon
,:line
or:point
geometries. using points or lines as polygons may have unexpected results.boundary
: for polygons, include pixels where the:center
is inside the polygon, where the polygon:touches
the pixel, or that are completely:inside
the polygon. The default is:center
.
fill
: the value or values to fill a polygon with. ASymbol
or tuple ofSymbol
will be used to retrieve properties from features or column values from table rows. An array or other iterable will be used for each geometry, in order.fill
can also be a function of the current value, e.g.x -> x + 1
.op
: A reducing function that accepts two values and returns one, likemin
tominimum
. For common methods this will be assigned for you, or is not required. But you can use it instead of areducer
as it will usually be faster.shape
: forcedata
to be treated as:polygon
,:line
or:point
, where possible Points can't be treated as lines or polygons, and lines may not work as polygons, but an attempt will be made.geometrycolumn
:Symbol
to manually select the column the geometries are in whendata
is a Tables.jl compatible table, or a tuple ofSymbol
for columns of point coordinates.progress
: show a progress bar,true
by default,false
to hide..verbose
: print information and warnings whne there are problems with the rasterisation.true
by default.threaded
: run operations in parallel.true
by default.filename
: a filename to write to directly, useful for large files.suffix
: a string or value to append to the filename. A tuple ofsuffix
will be applied to stack layers.keys(stack)
are the default.
Example
Rasterize a shapefile for China and plot, with a border.
using Rasters, RasterDataSources, ArchGDAL, Plots, Dates, Shapefile, Downloads
+using Rasters.LookupArrays
+
+# Download a borders shapefile
+shapefile_url = "https://github.com/nvkelso/natural-earth-vector/raw/master/10m_cultural/ne_10m_admin_0_countries.shp"
+shapefile_name = "country_borders.shp"
+isfile(shapefile_name) || Downloads.download(shapefile_url, shapefile_name)
+
+# Load the shapes for china
+china_border = Shapefile.Handle(shapefile_name).shapes[10]
+
+# Rasterize the border polygon
+china = rasterize(last, china_border; res=0.1, missingval=0, fill=1, boundary=:touches, progress=false)
+
+# And plot
+p = plot(china; color=:spring, legend=false)
+plot!(p, china_border; fillalpha=0, linewidth=0.6)
+
+savefig("build/china_rasterized.png"); nothing
+
+# output
+
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.rasterize!
— Functionrasterize!([reducer], dest, data; kw...)
Rasterize the geometries in data
into the Raster
or RasterStack
dest
, using the values specified by fill
.
Arguments
dest
: aRaster
orRasterStack
to rasterize into.reducer
: a reducing function to reduce the fill value for all geometries that cover or touch a pixel down to a single value. The default islast
. Any that takes an iterable and returns a single value will work, including custom functions. However, there are optimisations for built-in methods includingsum
,first
,last
,minimum
,maximum
,extrema
andStatistics.mean
. These may be an order of magnitude or more faster thancount
is a special-cased as it does not need a fill value.data
: a GeoInterface.jlAbstractGeometry
, or a nestedVector
ofAbstractGeometry
, or a Tables.jl compatible object containing a:geometry
column or points and values columns.
Keywords
These are detected automatically from A
and data
where possible.
fill
: the value or values to fill a polygon with. ASymbol
or tuple ofSymbol
will be used to retrieve properties from features or column values from table rows. An array or other iterable will be used for each geometry, in order.fill
can also be a function of the current value, e.g.x -> x + 1
.op
: A reducing function that accepts two values and returns one, likemin
tominimum
. For common methods this will be assigned for you, or is not required. But you can use it instead of areducer
as it will usually be faster.shape
: forcedata
to be treated as:polygon
,:line
or:point
, where possible Points can't be treated as lines or polygons, and lines may not work as polygons, but an attempt will be made.geometrycolumn
:Symbol
to manually select the column the geometries are in whendata
is a Tables.jl compatible table, or a tuple ofSymbol
for columns of point coordinates.progress
: show a progress bar,true
by default,false
to hide..verbose
: print information and warnings whne there are problems with the rasterisation.true
by default.threaded
: run operations in parallel.true
by default.to
: aRaster
,RasterStack
,Tuple
ofDimension
orExtents.Extent
. If noto
object is provided the extent will be calculated from the geometries, Additionally, when noto
object or anExtent
is passed forto
, thesize
orres
keyword must also be used.res
: the resolution of the dimensions, aReal
orTuple{<:Real,<:Real}
. Only required whento
is not used or is anExtents.Extent
, andsize
is not used.size
: the size of the output array, as aTuple{Int,Int}
or singleInt
for a square. Only required whento
is not used or is anExtents.Extent
, andres
is not used.crs
: acrs
which will be attached to the resulting raster whento
not passed or is anExtent
. Otherwise the crs fromto
is used.shape
: Forcedata
to be treated as:polygon
,:line
or:point
geometries. using points or lines as polygons may have unexpected results.boundary
: for polygons, include pixels where the:center
is inside the polygon, where the polygon:touches
the pixel, or that are completely:inside
the polygon. The default is:center
.
Example
using Rasters, RasterDataSources, ArchGDAL, Plots, Dates, Shapefile, GeoInterface, Downloads
+using Rasters.LookupArrays
+
+# Download a borders shapefile
+shapefile_url = "https://github.com/nvkelso/natural-earth-vector/raw/master/10m_cultural/ne_10m_admin_0_countries.shp"
+shapefile_name = "country_borders.shp"
+isfile(shapefile_name) || Downloads.download(shapefile_url, shapefile_name)
+
+# Load the shapes for indonesia
+indonesia_border = Shapefile.Handle(shapefile_name).shapes[1]
+
+# Make an empty EPSG 4326 projected Raster of the area of Indonesia
+dimz = X(Projected(90.0:0.1:145; sampling=Intervals(), crs=EPSG(4326))),
+ Y(Projected(-15.0:0.1:10.9; sampling=Intervals(), crs=EPSG(4326)))
+
+A = zeros(UInt32, dimz; missingval=UInt32(0))
+
+# Rasterize each indonesian island with a different number. The islands are
+# rings of a multi-polygon, so we use `GI.getring` to get them all separately.
+islands = collect(GeoInterface.getring(indonesia_border))
+rasterize!(last, A, islands; fill=1:length(islands), progress=false)
+
+# And plot
+p = plot(Rasters.trim(A); color=:spring)
+plot!(p, indonesia_border; fillalpha=0, linewidth=0.7)
+
+savefig("build/indonesia_rasterized.png"); nothing
+
+# output
+
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.resample
— Functionresample(x; kw...)
+resample(xs...; to=first(xs), kw...)
resample
uses ArchGDAL.gdalwarp
to resample a Raster
or RasterStack
to a new resolution
and optionally new crs
, or to snap to the bounds, resolution and crs of the object to
.
Arguments
x
: the object/s to resample.
Keywords
to
: aRaster
,RasterStack
,Tuple
ofDimension
orExtents.Extent
. If noto
object is provided the extent will be calculated fromx
,res
: the resolution of the dimensions, aReal
orTuple{<:Real,<:Real}
. Only required whento
is not used or is anExtents.Extent
, andsize
is not used.size
: the size of the output array, as aTuple{Int,Int}
or singleInt
for a square. Only required whento
is not used or is anExtents.Extent
, andres
is not used.crs
: acrs
which will be attached to the resulting raster whento
not passed or is anExtent
. Otherwise the crs fromto
is used.method
: ASymbol
orString
specifying the method to use for resampling. From the docs forgdalwarp
::near
: nearest neighbour resampling (default, fastest algorithm, worst interpolation quality).:bilinear
: bilinear resampling.:cubic
: cubic resampling.:cubicspline
: cubic spline resampling.:lanczos
: Lanczos windowed sinc resampling.:average
: average resampling, computes the weighted average of all non-NODATA contributing pixels. rms root mean square / quadratic mean of all non-NODATA contributing pixels (GDAL >= 3.3):mode
: mode resampling, selects the value which appears most often of all the sampled points.:max
: maximum resampling, selects the maximum value from all non-NODATA contributing pixels.:min
: minimum resampling, selects the minimum value from all non-NODATA contributing pixels.:med
: median resampling, selects the median value of all non-NODATA contributing pixels.:q1
: first quartile resampling, selects the first quartile value of all non-NODATA contributing pixels.:q3
: third quartile resampling, selects the third quartile value of all non-NODATA contributing pixels.:sum
: compute the weighted sum of all non-NODATA contributing pixels (since GDAL 3.1)
Where NODATA values are set to
missingval
.filename
: a filename to write to directly, useful for large files.suffix
: a string or value to append to the filename. A tuple ofsuffix
will be applied to stack layers.keys(stack)
are the default.
Note:
- GDAL may cause some unexpected changes in the data, such as returning a reversed Y dimension or changing the
crs
type fromEPSG
toWellKnownText
(it will represent the same CRS).
Example
Resample a WorldClim layer to match an EarthEnv layer:
using Rasters, RasterDataSources, ArchGDAL, Plots
+A = Raster(WorldClim{Climate}, :prec; month=1)
+B = Raster(EarthEnv{HabitatHeterogeneity}, :evenness)
+
+a = plot(A)
+b = plot(resample(A; to=B))
+
+savefig(a, "build/resample_example_before.png");
+savefig(b, "build/resample_example_after.png"); nothing
+# output
Before resample
:
After resample
:
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.replace_missing
— Functionreplace_missing(a::AbstractRaster, newmissingval)
+replace_missing(a::AbstractRasterStack, newmissingval)
Replace missing values in the array or stack with a new missing value, also updating the missingval
field/s.
Keywords
filename
: a filename to write to directly, useful for large files.suffix
: a string or value to append to the filename. A tuple ofsuffix
will be applied to stack layers.keys(stack)
are the default.
Example
using Rasters, RasterDataSources, ArchGDAL
+A = Raster(WorldClim{Climate}, :prec; month=1) |> replace_missing
+missingval(A)
+# output
+missing
Rasters.reproject
— Functionreproject(target::GeoFormat, x)
Reproject the dimensions of x
to a different crs.
Arguments
target
: any crs in a GeoFormatTypes.jl wrapper, e.g.EPSG
,WellKnownText
,ProjString
.x
: aDimension
,Tuple
ofDimension
,Raster
orRasterStack
.
Dimensions without an AbstractProjected
lookup (such as a Ti
dimension) are silently returned without modification.
reproject(source::GeoFormat, target::GeoFormat, dim::Dimension, val)
reproject
uses ArchGDAL.reproject, but implemented for a reprojecting a value array of values, a single dimension at a time.
Rasters.setcrs
— Functionsetcrs(x, crs)
Set the crs of a Raster
, RasterStack
, Tuple
of Dimension
, or a Dimension
. The crs
is expected to be a GeoFormatTypes.jl CRS
or Mixed
GeoFormat
type
Rasters.setmappedcrs
— Functionsetmappedcrs(x, crs)
Set the mapped crs of a Raster
, a RasterStack
, a Tuple
of Dimension
, or a Dimension
. The crs
is expected to be a GeoFormatTypes.jl CRS
or Mixed
GeoFormat
type
Base.skipmissing
— Functionskipmissing(itr::Raster)
Returns an iterable over the elements in a Raster
object, skipping any values equal to either the missingval
or missing
.
Rasters.trim
— Functiontrim(x; dims::Tuple, pad::Int)
Trim missingval(x)
from x
for axes in dims
, returning a view of x
.
Arguments
x
: ARaster
orRasterStack
. For stacks, all layers must having missing values for a pixel for it to be trimmed.
Keywords
dims
: By defaultdims=(XDim, YDim)
, so that trimming keeps the area ofX
andY
that contains non-missing values along all other dimensions.pad
: The trimmed size will be padded bypad
on all sides, although padding will not be added beyond the original extent of the array.
As trim
is lazy, filename
and suffix
keywords are not used.
Example
Create trimmed layers of Australian habitat heterogeneity.
using Rasters, RasterDataSources, Plots
+layers = (:evenness, :range, :contrast, :correlation)
+st = RasterStack(EarthEnv{HabitatHeterogeneity}, layers)
+
+# Roughly cut out australia
+ausbounds = X(100 .. 160), Y(-50 .. -10)
+aus = st[ausbounds...]
+a = plot(aus)
+
+# Trim missing values and plot
+b = plot(trim(aus))
+
+savefig(a, "build/trim_example_before.png");
+savefig(b, "build/trim_example_after.png"); nothing
+
+# output
+
Before trim
:
After trim
:
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.warp
— Functionwarp(A::AbstractRaster, flags::Dict; kw...)
Gives access to the GDALs gdalwarp
method given a Dict
of flag => value
arguments that can be converted to strings, or vectors where multiple space-separated arguments are required.
Arrays with additional dimensions not handled by GDAL (other than X
, Y
, Band
) are sliced, warped, and then combined to match the original array dimensions. These slices will not be written to disk and loaded lazyily at this stage - you will need to do that manually if required.
See the gdalwarp docs for a list of arguments.
Keywords
filename
: a filename to write to directly, useful for large files.suffix
: a string or value to append to the filename. A tuple ofsuffix
will be applied to stack layers.keys(stack)
are the default.
Any additional keywords are passed to ArchGDAL.Dataset
.
Example
This simply resamples the array with the :tr
(output file resolution) and :r
flags, giving us a pixelated version:
using Rasters, RasterDataSources, Plots
+A = Raster(WorldClim{Climate}, :prec; month=1)
+a = plot(A)
+
+flags = Dict(
+ :tr => [2.0, 2.0],
+ :r => :near,
+)
+b = plot(warp(A, flags))
+
+savefig(a, "build/warp_example_before.png");
+savefig(b, "build/warp_example_after.png"); nothing
+
+# output
+
Before warp
:
After warp
:
In practise, prefer resample
for this. But warp
may be more flexible.
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.zonal
— Functionzonal(f, x::Union{Raster,RasterStack}; of, kw...)
Calculate zonal statistics for the the zone of a Raster
or RasterStack
covered by the of
object/s.
Arguments
f
: any function that reduces an iterable to a single value, such assum
orStatistics.mean
x
: ARaster
orRasterStack
of
: ARaster
,RasterStack
, dim tuple, extent, GeoInterface.jl compatible geometry, Tables.jl compatible table of a:geometry
column, or anAbstractVector
of any of these objects..
Keywords
These can be used when of
is a GeoInterface.jl compatible object:
shape
: Forcedata
to be treated as:polygon
,:line
or:point
, where possible.boundary
: for polygons, include pixels where the:center
is inside the polygon, where the line:touches
the pixel, or that are completely:inside
inside the polygon. The default is:center
.progress
: show a progress bar,true
by default,false
to hide..
Example
``jldoctest using Rasters, RasterDataSources, ArchGDAL, Shapefile, DataFrames, Downloads, Statistics, Dates
Download a borders shapefile
neurl = "https://github.com/nvkelso/natural-earth-vector/raw/master/10mcultural/ne10madmin0countries" shpurl, dbfurl = neurl * ".shp", neurl * ".dbf" shpname, dbfname = "countryborders.shp", "countryborders.dbf" isfile(shpname) || Downloads.download(shpurl, shpname) isfile(dbfurl) || Downloads.download(dbfurl, dbfname)
Download and read a raster stack from WorldClim
st = RasterStack(WorldClim{Climate}; month=Jan, lazy=false)
Load the shapes for world countries
countries = Shapefile.Table(shp_name) |> DataFrame
Calculate the january mean of all climate variables for all countries
january_stats = zonal(mean, st; of=countries, boundary=:touches, progress=false) |> DataFrame
Add the country name column (natural earth has some string errors it seems)
insertcols!(january_stats, 1, :country => first.(split.(countries.ADMIN, r"[^A-Za-z ]")))
output
258×8 DataFrame Row │ country tmin tmax tavg prec ⋯ │ SubStrin… Float32 Float32 Float32 Float64 ⋯ ─────┼────────────────────────────────────────────────────────────────────────── 1 │ Indonesia 21.5447 29.1864 25.3656 271.063 ⋯ 2 │ Malaysia 21.3087 28.4291 24.8688 273.381 3 │ Chile 7.24534 17.9263 12.5858 78.1287 4 │ Bolivia 17.2065 27.7454 22.4759 192.542 5 │ Peru 15.0273 25.5504 20.2888 180.007 ⋯ 6 │ Argentina 13.6751 27.6715 20.6732 67.1837 7 │ Dhekelia Sovereign Base Area 5.87126 15.8991 10.8868 76.25 8 │ Cyprus 5.65921 14.6665 10.1622 97.4474 ⋮ │ ⋮ ⋮ ⋮ ⋮ ⋮ ⋱ 252 │ Spratly Islands 25.0 29.2 27.05 70.5 ⋯ 253 │ Clipperton Island 21.5 33.2727 27.4 6.0 254 │ Macao S 11.6694 17.7288 14.6988 28.0 255 │ Ashmore and Cartier Islands NaN NaN NaN NaN 256 │ Bajo Nuevo Bank NaN NaN NaN NaN ⋯ 257 │ Serranilla Bank NaN NaN NaN NaN 258 │ Scarborough Reef NaN NaN NaN NaN 3 columns and 243 rows omitted ```
Slice and combine to and from RasterSeries
Rasters.slice
— Functionslice(A::Union{AbstractRaster,AbstractRasterStack,AbstracRasterSeries}, dims) => RasterSeries
Slice views along some dimension/s to obtain a RasterSeries
of the slices.
For a Raster
or RasterStack
this will return a RasterSeries
of Raster
or RasterStack
that are slices along the specified dimensions.
For a RasterSeries
, the output is another series where the child objects are sliced and the series dimensions index is now of the child dimensions combined. slice
on a RasterSeries
with no dimensions will slice along the dimensions shared by both the series and child object.
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
Rasters.combine
— Functioncombine(A::Union{AbstractRaster,AbstractRasterStack,AbstracRasterSeries}, [dims]) => Raster
Combine a RasterSeries
along some dimension/s, creating a new Raster
or RasterStack
, depending on the contents of the series.
If dims
are passed, only the specified dimensions will be combined with a RasterSeries
returned, unless dims
is all the dims in the series.
WARNING: This feature is experimental. It may change in future versions, and may not be 100% reliable in all cases. Please file github issues if problems occur.
File operations
These Base
and DimensionalData
methods have specific Rasters.jl versions:
DimensionalData.modify
— Functionmodify(f, series::AbstractRasterSeries)
Apply function f
to the data of the child object. If the child is an AbstractRasterStack
the function will be passed on to its child AbstractRaster
s.
f
must return an idenically sized array.
This method triggers a complete rebuild of all objects, and disk based objects will be transferred to memory.
An example of the usefulnesss of this is for swapping out array backend for an entire series to CuArray
from CUDA.jl to copy data to a GPU.
Base.open
— Functionopen(f, A::AbstractRaster; write=false)
open
is used to open any lazy=true
AbstractRaster
and do multiple operations on it in a safe way. The write
keyword opens the file in write lookup so that it can be altered on disk using e.g. a broadcast.
f
is a method that accepts a single argument - an Raster
object which is just an AbstractRaster
that holds an open disk-based object. Often it will be a do
block:
lazy=false
(in-memory) rasters will ignore open
and pass themselves to f
.
# A is an `Raster` wrapping the opened disk-based object.
+open(Raster(filepath); write=true) do A
+ mask!(A; with=maskfile)
+ A[I...] .*= 2
+ # ... other things you need to do with the open file
+end
By using a do block to open files we ensure they are always closed again after we finish working with them.
Base.read
— Functionread(A::AbstractRaster)
+read(A::AbstractRasterStack)
+read(A::AbstractRasterSeries)
read
will move a Rasters.jl object completely to memory.
Base.read!
— Functionread!(src::Union{AbstractString,AbstractRaster}, dst::AbstractRaster)
+read!(src::Union{AbstractString,AbstractRasterStack}, dst::AbstractRasterStack)
+read!(scr::AbstractRasterSeries, dst::AbstractRasterSeries)
read!
will copy the data from src
to the object dst
.
src
can be an object or a file-path String
.
Base.write
— FunctionBase.write(filename::AbstractString, A::AbstractRaster; kw...)
Write an AbstractRaster
to file, guessing the backend from the file extension.
Keyword arguments are passed to the write
method for the backend.
Base.write(filename::AbstractString, s::AbstractRasterStack; suffix, kw...)
Write any AbstractRasterStack
to file, guessing the backend from the file extension.
Keywords
suffix
: suffix to append to file names. By default the layer key is used.
Other keyword arguments are passed to the write
method for the backend.
If the source can't be saved as a stack-like object, individual array layers will be saved.
Base.write(filepath::AbstractString, s::AbstractRasterSeries; kw...)
Write any AbstractRasterSeries
to file, guessing the backend from the file extension.
The lookup values of the series will be appended to the filepath (before the extension), separated by underscores.
Keywords
See other docs for write
. All keywords are passed through to Raster
and RasterStack
methods.
Base.write(filename::AbstractString, ::Type{GRDsource}, s::AbstractRaster; force=false)
Write a Raster
to a .grd file with a .gri header file. The extension of filename
will be ignored.
Returns filename
.
Base.write(filename::AbstractString, ::Type{NCDsource}, A::AbstractRaster)
Write an NCDarray to a NetCDF file using NCDatasets.jl
Returns filename
.
Base.write(filename::AbstractString, ::Type{NCDsource}, s::AbstractRasterStack; kw...)
Write an NCDstack to a single netcdf file, using NCDatasets.jl.
Currently Metadata
is not handled for dimensions, and Metadata
from other AbstractRaster
@types is ignored.
Keywords
Keywords are passed to NCDatasets.defVar
.
append
: If true, the variable of the current Raster will be appended tofilename
. Note that the variable of the current Raster should be not exist before. If not, you need to setappend = false
. Rasters.jl can not overwrite a previous existing variable.fillvalue
: A value filled in the NetCDF file to indicate missing data. It will be stored in the_FillValue
attribute.chunksizes
: Vector integers setting the chunk size. The total size of a chunk must be less than 4 GiB.deflatelevel
: Compression level: 0 (default) means no compression and 9 means maximum compression. Each chunk will be compressed individually.shuffle
: If true, the shuffle filter is activated which can improve the compression ratio.checksum
: The checksum method can be:fletcher32
or:nochecksum
(checksumming is disabled, which is the default)typename
(string): The name of the NetCDF type required for vlen arrays (https://web.archive.org/save/https://www.unidata.ucar.edu/software/netcdf/netcdf-4/newdocs/netcdf-c/nc005fdef005fvlen.html)
Base.write(filename::AbstractString, ::Type{GDALsource}, A::AbstractRaster; force=false, kw...)
Write a Raster
to file using GDAL.
Keywords
driver
: A GDAL driver name or a GDAL driver retrieved viaArchGDAL.getdriver(drivername)
. Guessed from the filename extension by default.options::Dict{String,String}
: A dictionary containing the dataset creation options passed to the driver. For example:Dict("COMPRESS"=>"DEFLATE")
Valid options for the drivers can be looked up here: https://gdal.org/drivers/raster/index.html
Returns filename
.
Internals
Rasters.FileArray
— TypeFileArray{X} <: DiskArrays.AbstractDiskArray
Filearray is a DiskArrays.jl AbstractDiskArray
. Instead of holding an open object, it just holds a filename string that is opened lazily when it needs to be read.
Rasters.FileStack
— TypeFileStack{X,K}
+
+FileStack{X,K}(filename, types, sizes, eachchunk, haschunks, write)
A wrapper object that holds file pointer and size/chunking metadata for a multi-layered stack stored in a single file, typically netcdf or hdf5.
X
is a backend type like NCDsource
, and K
is a tuple of Symbol
keys.
Rasters.OpenStack
— TypeOpenStack{X,K}
+
+OpenStack{X,K}(dataset)
A wrapper for any stack-like opened dataset that can be indexed with Symbol
keys to retrieve AbstractArray
layers.
OpenStack
is usually hidden from users, wrapped in a regular RasterStack
passed as the function argument in open(stack)
when the stack is contained in a single file.
X
is a backend type like NCDsource
, and K
is a tuple of Symbol
keys.
Rasters.RasterDiskArray
— TypeRasterDiskArray <: DiskArrays.AbstractDiskArray
A basic DiskArrays.jl wrapper for objects that don't have one defined yet. When we open
a FileArray
it is replaced with a RasterDiskArray
.