IceFloeTracker.jl
Overview
IceFloeTracker.jl is a collection of routines and tools for processing remote sensing imagery, identifying sea ice floes, and tracking the displacement and rotation of ice floes across multiple images. It can be used either standalone to create custom processing pathways or with the Ice Floe Tracker Pipeline.
Algorithm components
The Ice Floe Tracker (IFT) package includes the core functions for the three main steps of the algorithm. These functions can be used independently and can be customized for specific use cases.
Preprocessing
IFT operates on optical satellite imagery. The main functions are designed with "true color" and "false color" imagery in mind, and have thus far primarily been tested on imagery from the Moderate Resolution Imaging Spectroradiometer (MODIS) from the NASA Aqua and Terra satellites. The preprocessing routines mask land and cloud features, and aim to adjust and sharpen the remainder of the images to amplify the contrast along the edges of sea ice floes. (TBD: Link to main preprocessing page)
Segmentation
The IFT segmentation functions include functions for semantic segmentation (pixel-by-pixel assignment into predefined categories) and object-based segmentation (groupings of pixels into distinct objects). The semantic segmentation steps use $k$-means to group pixels into water and ice regions. A combination of watershed functions, morphological operations, and further applications of $k$-means are used to identify candidate ice floes. (TBD: Link to main segmentation page)
Tracking
Ice floe tracking is carried out by comparing the shapes produced in the segmentation step. Shapes with similar area are rotated until the difference in surface area is minimized, and then the edge shapes are compared using a Ѱ-s curve. If thresholds for correlation and area differences are met, then the floe with the best correlation and smallest area differences are considered matches and the objects are assigned the same label. In the end, trajectories for individual floes are recorded in a dataframe.
Developers
IceFloeTracker.jl is a product of the Wilhelmus Lab at Brown University, led by Monica M. Wilhelmus. The original algorithm was developed by Rosalinda Lopez-Acosta during her PhD work at University of California Riverside, advised by Dr. Wilhelmus. The translation of the original Matlab code into the current modular, open source Julia package has been carried out in conjunction with the Center for Computing and Visualization at Brown University. Contributors include Daniel Watkins, Maria Isabel Restrepo, Carlos Paniagua, Tim Divoll, John Holland, and Bradford Roarr.
Citing
If you use IceFloeTracker.jl in research, teaching, or elsewhere, please mention the IceFloeTracker package and cite our journal article outlining the algorithm:
Lopez-Acosta et al., (2019). Ice Floe Tracker: An algorithm to automatically retrieve Lagrangian trajectories via feature matching from moderate-resolution visual imagery. Remote Sensing of Environment, 234(111406), doi:10.1016/j.rse.2019.111406.
Papers using Ice Floe Tracker
- Manucharyan, Lopez-Acosta, and Wilhelmus (2022)*. Spinning ice floes reveal intensification of mesoscale eddies in the western Arctic Ocean. Scientific Reports, 12(7070), doi:10.1038/s41598-022-10712-z
- Watkins, Bliss, Hutchings, and Wilhelmus (2023)*. Evidence of Abrupt Transitions Between Sea Ice Dynamical Regimes in the East Greenland Marginal Ice Zone. Geophysical Research Letters, 50(e2023GL103558), pp. 1-10, doi:10.1029/2023GL103558
*Papers using data from the Matlab implementation of Ice Floe Tracker.
Functions
IceFloeTracker.add_passtimes!
— Methodadd_passtimes!(props, passtimes)
Add a column passtime
to each DataFrame in props
containing the time of the image in which the floes were captured.
Arguments
props
: array of DataFrames containing floe properties.passtimes
: array ofDateTime
objects containing the time of the image in which the floes were captured.
IceFloeTracker.addfloemasks!
— Methodaddfloemasks!(props::DataFrame, floeimg::FloeLabelsImage)
Add a column to props
called floearray
containing the cropped floe masks from floeimg
.
IceFloeTracker.addlatlon!
— Methodaddlatlon(pairedfloesdf::DataFrame, refimage::AbstractString)
Add columns latitude
, longitude
, and pixel coordinates x
, y
to pairedfloesdf
.
Arguments
pairedfloesdf
: dataframe containing floe tracking data.refimage
: path to reference image.
IceFloeTracker.addψs!
— Methodaddψs!(props::Vector{DataFrame})
Add the ψ-s curves to each member of props
.
Note: each member of props
must have a mask
column with a binary image representing the floe.
To add floe masks see addfloemasks!
.
IceFloeTracker.branch
— Methodbranch(img::AbstractArray{Bool})
Find branch points in skeletonized image img
according to Definition 3 of [1].
[1] Arcelli, Carlo, and Gabriella Sanniti di Baja. "Skeletons of planar patterns." Machine Intelligence and Pattern Recognition. Vol. 19. North-Holland, 1996. 99-143.
IceFloeTracker.bridge
— Methodbridge(bw)
Set 0-valued pixels to 1 if they have two nonzero neighbors that are not connected. Note the following exceptions:
0 0 0 0 0 0 1 0 1 becomes 1 1 1 0 0 0 0 0 0
1 0 1 1 1 1 0 0 0 becomes 0 0 0 0 0 0 0 0 0
The same applies to all their corresponding rotations.
Examples
julia> bw = [0 0 0; 0 0 0; 1 0 1]
3×3 Matrix{Int64}:
0 0 0
0 0 0
1 0 1
julia> bridge(bw)
3×3 BitMatrix:
0 0 0
0 0 0
1 1 1
julia> bw = [1 0 0; 1 0 1; 0 0 1]
3×3 Matrix{Int64}:
1 0 0
1 0 1
0 0 1
julia> bridge(bw)
3×3 BitMatrix:
1 1 0
1 1 1
0 1 1
julia> bw = [1 0 1; 0 0 0; 1 0 1]
3×3 Matrix{Int64}:
1 0 1
0 0 0
1 0 1
julia> bridge(bw)
3×3 BitMatrix:
1 1 1
1 1 1
1 1 1
IceFloeTracker.convertcentroid!
— Methodconvertcentroid!(propdf, latlondata, colstodrop)
Convert the centroid coordinates from row and column to latitude and longitude dropping unwanted columns specified in colstodrop
for the output data structure. Addionally, add columns x
and y
with the pixel coordinates of the centroid.
IceFloeTracker.converttounits!
— Methodconverttounits!(propdf, latlondata, colstodrop)
Convert the floe properties from pixels to kilometers and square kilometers where appropiate. Also drop the columns specified in colstodrop
.
IceFloeTracker.create_cloudmask
— Functioncreate_cloudmask(img; f=AbstractCloudMaskAlgorithm)
Cloud masks in the IFT are BitMatrix objects such that for an image I and cloudmask C, cloudy pixels can be selected by I[C], and clear-sky pixels can be selected with I[.!C]. Construction of a cloud mask uses the syntax
f = CloudMaskAlgorithm(parameters)
C = create_cloudmask(img; CloudMaskAlgorithm)
By default, create_cloudmask
uses the algorithm found in [1]. This algorithm converts a 3-channel MODIS 7-2-1 false color image into a 1-channel binary matrix in which clouds = 1 and anything else = 0. The algorithm aims to identify patches of opaque cloud while allowing thin and transparent cloud to remain. This algorithm is instantiated using
f = LopezAcostaCloudMask()
In this case, the default values are applied. It can also called using a set of customized parameters. These values must be real numbers between 0 and 1. To reproduce the default parameters, you may call
f = LopezAcostaCloudMask(prelim_threshold=110/255, band_7_threshold=200/255, band_2_threshold=190/255, ratio_lower=0.0, ratio_upper=0.75).
A stricter cloud mask was defined in [2], covering more cloudy pixels while minimally impacting the masking of cloud-covered ice pixels.
f = LopezAcostaCloudMask(prelim_threshold=53/255, band_7_threshold=130/255, band_2_threshold=169/255, ratio_lower=0.0, ratio_upper=0.53).
These parameters together define a piecewise linear partition of pixels based on their Band 7 and Band 2 callibrated reflectance. Pixels with intensity above prelim_threshold
are considered as potential cloudy pixels. Then, pixels with Band 7 reflectance less than band_7_threshold
, Band 2 reflectance greater than band_2_threshold
, and Band 7 to Band 2 ratios between ratio_lower
and ratio_upper
are removed from the cloud mask (i.e., set to cloud-free).
Arguments
false_color_image
: corrected reflectance false color image - bands [7,2,1]prelim_threshold
: threshold value used to identify clouds in band 7, N0f8(RGB intensity/255)band_7_threshold
: threshold value used to identify cloud-ice in band 7, N0f8(RGB intensity/255)band_2_threshold
: threshold value used to identify cloud-ice in band 2, N0f8(RGB intensity/255)ratio_lower
: threshold value used to set lower ratio of cloud-ice in bands 7 and 2ratio_upper
: threshold value used to set upper ratio of cloud-ice in bands 7 and 2ratio_offset
: offset value used to adjust the upper ratio of cloud-ice in bands 7 and 2
- Lopez-Acosta, R., Schodlok, M. P., & Wilhelmus, M. M. (2019). Ice Floe Tracker: An algorithm to automatically retrieve Lagrangian trajectories via feature matching from moderate-resolution visual imagery. Remote Sensing of Environment, 234(111406), 1–15. (https://doi.org/10.1016/j.rse.2019.111406)[https://doi.org/10.1016/j.rse.2019.111406]
- Watkins, D.M., Kim, M., Paniagua, C., Divoll, T., Holland, J.G., Hatcher, S., Hutchings, J.K., and Wilhelmus, M.M. (in prep). Calibration and validation of the Ice Floe Tracker algorithm.
IceFloeTracker.create_landmask
— Methodcreate_landmask(landmask_image, struct_elem, fill_value_lower, fill_value_upper)
Convert a 3-channel RGB land mask image to a 1-channel binary matrix, and use a structuring element to extend a buffer to mask complex coastal features. In the resulting mask, land = 0 and ocean = 1. Returns a named tuple with the dilated and non-dilated landmask.
Arguments
landmask_image
: RGB land mask image fromfetchdata
struct_elem
: structuring element for dilation (optional)fill_value_lower
: fill holes having at least these many pixels (optional)fill_value_upper
: fill holes having at most these many pixels (optional)
IceFloeTracker.cropfloe
— Methodcropfloe(floesimg, props, i)
Crops the floe delimited by the bounding box data in props
at index i
from the floe image floesimg
.
If the dataframe has bounding box data min_row
, min_col
, max_row
, max_col
, but no label
, then returns the largest contiguous component.
If the dataframe has bounding box data min_row
, min_col
, max_row
, max_col
, and a label
, then returns the component with the label. In this case, floesimg
must be an Array{Int}.
If the dataframe has only a label
and no bounding box data, then returns the component with the label, padded by one cell of zeroes on all sides. In this case, floesimg
must be an Array{Int}.
IceFloeTracker.cropfloe
— Methodcropfloe(floesimg, min_row, min_col, max_row, max_col)
Crops the floe delimited by min_row
, min_col
, max_row
, max_col
, from the floe image floesimg
.
IceFloeTracker.cropfloe
— Methodcropfloe(floesimg, min_row, min_col, max_row, max_col, label)
Crops the floe from floesimg
with the label label
, returning the region bounded by min_row
, min_col
, max_row
, max_col
, and converting to a BitMatrix.
IceFloeTracker.get_rotation_measurements
— MethodCalculate the angle and rotation rate between two observations in DataFrameRows row1
and row2
. image_column
and time_column
specify which columns to use from the DataFrameRows. registration_function
is used to compare the two images and should return an angle. Returns a NamedTuple with the angle theta_rad
, time difference dt_sec
and rotation rate omega_rad_per_sec
, and the two input rows.
IceFloeTracker.get_rotation_measurements
— MethodCalculate the angle and rotation rate between a measurement in a DataFrameRow measurement
, and all the other rows in DataFrame df
.
image_column
is the column with the image to compare,time_column
is the column with the timepoint of each observation,registration_function
is used to compare the two images and should return an angle.
Returns a vector of NamedTuple
s with one entry for each comparison, with the angle theta_rad
, time difference dt_sec
and rotation rate omega_rad_per_sec
, and the two input rows for each comparison row1
and row2
.
IceFloeTracker.get_rotation_measurements
— MethodCalculate the angle and rotation rate between observations in DataFrame df
.
id_column
is the column with the ID of the image over several observations, e.g. the floe ID.image_column
is the column with the image to compare,time_column
is the column with the timepoint of each observation,registration_function
is used to compare the two images and should return an angle.
Each row is compared to each other row in df
which are:
- for the same object ID,
- strictly older,
- not older than the previous day.
Returns a DataFrame with one row for each comparison, with the angle theta_rad
, time difference dt_sec
and rotation rate omega_rad_per_sec
, and all the other values from df
with the column name suffix 1
for the first observation and 2
for the second.
IceFloeTracker.get_rotation_measurements
— MethodCalculate the angle and rotation rate between two images image1
and image2
at times time1
and time2
. Returns a NamedTuple with the angle theta_rad
, time difference dt_sec
and rotation rate omega_rad_per_sec
. registration_function
is used to compare the two images and should return an angle.
IceFloeTracker.imsharpen
— Functionimsharpen(truecolor_image, landmask_no_dilate, lambda, kappa, niters, nbins, rblocks, cblocks, clip, smoothing_param, intensity)
Sharpen truecolor_image
.
Arguments
truecolor_image
: input image in truecolorlandmask_no_dilate
: landmask for region of interestlambda
: speed of diffusion (0–0.25)kappa
: conduction coefficient for diffusion (25–100)niters
: number of iterations of diffusionnbins
: number of bins during histogram equalizationrblocks
: number of row blocks to divide input image during equalizationcblocks
: number of column blocks to divide input image during equalizationclip
: Thresholds for clipping histogram bins (0–1); values closer to one minimize contrast enhancement, values closer to zero maximize contrast enhancementsmoothing_param
: pixel radius for gaussian blurring (1–10)intensity
: amount of sharpening to perform
IceFloeTracker.latlon
— Methodlatlon(imgpath::AbstractString)
Reads the GeoTiff located at <imgpath>, extracts the coordinate reference system, and produces a lookup table with for the column and row values in the same projection as the GeoTiff and a 2D array for latitude and longitude.
IceFloeTracker.loadimg
— Methodloadimg(; dir::String, fname::String)
Load an image from dir
with filename fname
into a matrix of Float64
values. Returns the loaded image.
IceFloeTracker.long_tracker
— Methodlong_tracker(props, condition_thresholds, mc_thresholds)
Track ice floes over multiple observations.
Trajectories are built as follows:
- Assume the floes detected in observation 1 are trajectories of length 1.
- For each subsequent observation:
- Determine the latest observation for each trajectory – these are the "current trajectory heads".
- Find matches between the the current trajectory heads and the new observed floes, extending those trajectories.
- Any unmatched floe in an observation is added as a new trajectory starting point.
Arguments
props::Vector{DataFrame}
: A vector of DataFrames, each containing ice floe properties for a single day. Each DataFrame must have the following columns:- "area"
- "min_row"
- "min_col"
- "max_row"
- "max_col"
- "row_centroid"
- "col_centroid"
- "convex_area"
- "majoraxislength"
- "minoraxislength"
- "orientation"
- "perimeter"
- "mask": 2D array of booleans
- "passtime": A timestamp for the floe
- "psi": the psi-s curve for the floe
- "uuid": a universally unique identifier for each segmented floe
condition_thresholds
: namedtuple of thresholds for deciding whether to match floei
from dayk
to floe j from dayk+1
. SeeIceFloeTracker.condition_thresholds
for sample values.mc_thresholds
: thresholds for area mismatch and psi-s shape correlation. SeeIceFloeTracker.mc_thresholds
for sample values.
Returns
A DataFrame with the above columns, plus extra columns:
area_mismatch
andcorr
, which are the area mismatch and correlation between a floe and the one that preceeds it in the trajectory.head_uuid
, the floe which was best matched by this floe.- Trajectories are identified by:
- a unique identifier
ID
and the - UUID of the trajectory,
trajectory_uuid
.
- a unique identifier
IceFloeTracker.matchcorr
— Methodmatchcorr(
f1::T,
f2::T,
Δt::F,
mxrot::S=10,
psi::F=0.95,
sz::S=16,
comp::F=0.25,
mm::F=0.22
)
where {T<:AbstractArray{Bool,2},S<:Int64,F<:Float64}
Compute the mismatch mm
and psi-s-correlation c
for floes with masks f1
and f2
.
The criteria for floes to be considered equivalent is as follows: - c
greater than mm
- _mm
is less than mm
A pair of NaN
is returned for cases for which one of their mask dimension is too small or their sizes are not comparable.
Arguments
f1
: mask of floe 1f2
: mask of floe 2Δt
: time difference between floesmxrot
: maximum rotation (in degrees) allowed between floesn (default: 10)psi
: psi-s-correlation threshold (default: 0.95)sz
: size threshold (default: 16)comp
: size comparability threshold (default: 0.25)mm
: mismatch threshold (default: 0.22)
IceFloeTracker.padnhood
— Methodpadnhood(img, I, nhood)
Pad the matrix img[nhood]
with zeros according to the position of I
within the edgesimg
.
Returns img[nhood]
if I
is not an edge index.
IceFloeTracker.regionprops_table
— Functionregionprops_table(label_img, intensity_img; properties, connectivity, extra_properties)
A wrapper of the regionprops_table
function from the skimage python library.
See its full documentation at https://scikit-image.org/docs/stable/api/skimage.measure.html#regionprops-table.
Arguments
label_img
: Image with the labeled objects of interestintensity_img
: (Optional) Used for generatingextra_properties
, integer/float array from which (presumably)label_img
was generatedproperties
: List (Vector
orTuple
) of properties to be generated for each connected component inlabel_img
extra_properties
: (Optional) not yet implemented. It will be set tonothing
Notes
- Zero indexing has been corrected for the
bbox
andcentroid
properties bbox
data (max_col
andmax_row
) are inclusivecentroid
data are rounded to the nearest integer
See also regionprops
Examples
julia> using IceFloeTracker, Random
julia> Random.seed!(123);
julia> bw_img = rand([0, 1], 5, 10)
5×10 Matrix{Int64}:
1 0 1 0 0 0 0 0 0 1
1 0 1 1 1 0 0 0 1 1
1 1 0 1 1 0 1 0 0 1
0 1 0 1 0 0 0 0 1 0
1 0 0 0 0 1 0 1 0 1
julia> label_img = IceFloeTracker.label_components(bw_img, trues(3,3))
5×10 Matrix{Int64}:
1 0 1 0 0 0 0 0 0 4
1 0 1 1 1 0 0 0 4 4
1 1 0 1 1 0 3 0 0 4
0 1 0 1 0 0 0 0 4 0
1 0 0 0 0 2 0 4 0 4
julia> properties = ["area", "perimeter"]
2-element Vector{String}:
"area"
"perimeter"
julia> IceFloeTracker.regionprops_table(label_img, bw_img, properties = properties)
4×2 DataFrame
Row │ area perimeter
│ Int32 Float64
─────┼──────────────────
1 │ 13 11.6213
2 │ 1 0.0
3 │ 1 0.0
4 │ 7 4.62132
IceFloeTracker.@persist
— Macro@persist img fname
@persist(img,fname)
@persist img
@persist(img)
@persist img fname ts
@persist(img, fname, ts)
Given a reference to an image object img
, the macro persists (saves to a file) img
to the current working directory using fname
as filename. Returns img
.
Arguments
img
: Symbol expression representing an image object loaded in memory.fname
: Optional filename for the persisted image.ts
: Optional boolean to attach timestamp tofname
.
Index
IceFloeTracker.add_passtimes!
IceFloeTracker.addfloemasks!
IceFloeTracker.addlatlon!
IceFloeTracker.addψs!
IceFloeTracker.branch
IceFloeTracker.bridge
IceFloeTracker.convertcentroid!
IceFloeTracker.converttounits!
IceFloeTracker.create_cloudmask
IceFloeTracker.create_landmask
IceFloeTracker.cropfloe
IceFloeTracker.cropfloe
IceFloeTracker.cropfloe
IceFloeTracker.get_rotation_measurements
IceFloeTracker.get_rotation_measurements
IceFloeTracker.get_rotation_measurements
IceFloeTracker.get_rotation_measurements
IceFloeTracker.imsharpen
IceFloeTracker.latlon
IceFloeTracker.loadimg
IceFloeTracker.long_tracker
IceFloeTracker.matchcorr
IceFloeTracker.padnhood
IceFloeTracker.regionprops_table
IceFloeTracker.@persist