Skip to content

Simulation API

These are the full docstrings for the Simulation subsection of Sargassum.jl.

Index

Interpolants

Sargassum.GriddedField Type
julia
struct GriddedField{N, T, U}

A container for gridded data, possibly time-dependent.

Fields

  • dims_names: A Vector of Tuple{Symbol, Unitful.Unitlike}s such that dims_names[i][1] is the ith dimension of fields and dims_names[i][2] gives the units of the ith dimension.

  • dims: A Dict mapping variable names to ranges they take.

  • fields_names: A Vector of Tuple{Symbol, Unitful.Unitlike}s such that fields_names[i][1] is the ith field and fields_names[i][2] gives its units.

  • fields: A Dict mapping field names to their arrays.

Example

If GriddedField.dims_names == [(:x, u"km"), (:y, u"km"), (:t, u"d")], this implies that the order of the dimensions of each field is (x, y, t) and further that the units of x, y and t are km, km and d, respectively.

Constructor

GriddedField(n_dims; floats = Float64, ints = Int64)

where n_dims is the number of dimensions of the field and floats and ints give the datatypes used. For example a water velocity field would have n_dims = 3 from (x, y, t). A land interpolant would have n_dims = 2 from (x, y), i.e. the land location is not time-dependent.

source

Sargassum.InterpolatedField Type
julia
struct InterpolatedField{N, T, U, I}

A container for interpolants of gridded data, possibly time-dependent.

Fields

  • dims_names: A Vector of Tuple{Symbol, Unitful.Unitlike}s such that dims_names[i][1] is the ith dimension of fields and dims_names[i][2] gives the units of the ith dimension.

  • dims: A Dict mapping variable names to ranges they take.

  • fields_names: A Vector of Tuple{Symbol, Unitful.Unitlike}s such that fields_names[i][1] is the ith field and fields_names[i][2] gives its units.

  • fields: A Dict mapping field names to their interpolants.

Example

If InterpolatedField.dims_names == [(:x, u"km"), (:y, u"km"), (:t, u"d")], this implies that the order of the dimensions of each field is (x, y, t) and further that the units of x, y and t are km, km and d, respectively.

Constructor

InterpolatedField(gf; interpolant_type = "cubic", extrapolate_value = 0.0)

where gf is a GriddedField and with the optional arguments

  • interpolant_type: Two convenience flags are provided, "cubic" and "nearest" which refer to cubic BSpline and nearest-neighbor interpolation, respectively. Alternatively, any Interpolations.InterpolationType can be provided. Default "cubic".

  • extrapolate_value: A constant extrapolation is performed with this value. Default "0.0".

Plotting

This object can be viz and viz!.

source

Sargassum.add_derivatives! Method
julia
add_derivatives(itrf; interpolant_type, extrapolate_value, xyt_names, vxvy_names, Dx_Dy_vort_names, geometry)

Add three additional fields to InterpolatedField, namely the x and y components of the material derivative and the vorticity.

Arguments

  • itrf: An InterpolatedField. The field should have x, y and t variables along with x and y components of the corresponding field (e.g. water currents.)

Optional Arguments

  • interpolant_type: Two convenience flags are provided, "cubic" and "nearest" which refer to cubic BSpline and nearest-neighbor interpolation, respectively. Alternatively, any Interpolations.InterpolationType can be provided. Default "cubic".

  • extrapolate_value: A constant extrapolation is performed with this value. Default "0.0".

  • xyt_names: A Tuple with three symbols, corresponding to the x,y and t variables in that order. Default (:x, :y, :t).

  • vxvy_names: A Tuple with two symbols, corresponding to the x and y components of the vector field, in that order. Default (:u, :v).

  • Dx_Dy_vort_names: A Tuple with three symbols, corresponding to the x component of the material derivative, the y component of the material derivative and the vorticity, in that order. Default (:DDt_x, :DDt_y, :vorticity).

  • geometry: If true, include factors that take into account the spherical geometry of the Earth. Default true.

source

Sargassum.add_field! Method
julia
add_field!(gf, infile, field_name_in, field_name_out, field_units_in, field_units_out; take_axes, permutation, scale_factor_name, add_offset_name, missings_name, missings_replacement)

Add a new field to gf::GriddedField with data read from a NetCDF or MAT file infile.

The new dimension appears last in the list of field names.

Arguments

  • gf: The GriddedField to be modified.

  • infile: The path to the NetCDF/MAT file.

  • field_name_in: A String giving the name of the field to read in as it appears in the NetCDF/MAT file.

  • field_name_out: A Symbol giving the name of the added field in gf.

  • field_units_in: A Unitful.Unitlike giving the units of the field as they appear in the NetCDF/MAT file.

  • field_units_out: A String giving the kind of quantity being read; should be one of keys(UNITS).

Optional Arguments

  • take_axes: If provided, only the selected elements will be taken via field[take_axes...]. For example, if the field is four dimensional, passing take_axes = [:,:,1,:] would result in a three dimensional field with dimensions 1, 2 and 4 preserved - indexed on the first element of the third dimension. Default nothing.

  • permutation: If provided, the field will be permuted according to permutation. Applied AFTER take_axes. Default nothing.

  • scale_factor_name: The name of the scale factor, only for NetCDF files. If no scale factor is found, it is taken to be 1. Default "scale_factor".

  • add_offset_name: The name of the additive offset, only for NetCDF files. If no additive offset is found, it is taken to be 0. Default "add_offset".

  • missings_name: A vector of names of missing/fill/extra values, only for NetCDF files. Each such value will be replaced by missings_replacement if found. Default ["_FillValue", "missing_value"].

  • missings_replacement: missings_name replaces the missing/fill/extra values with this. Default 0.0.

source

Sargassum.add_spatial_dimension! Method
julia
add_spatial_dimension!(gf, infile, dim_name_in, dim_name_out, dim_units_in, dim_units_out; transform)

Add a new spatial dimension to gf::GriddedField with data read from a NetCDF or MAT file infile.

The new dimension appears last in the list of dimension names.

Arguments

  • gf: The GriddedField to be modified.

  • infile: The path to the NetCDF/MAT file.

  • dim_name_in: A String giving the name of the dimension to read in as it appears in the NetCDF/MAT file.

  • dim_name_out: A Symbol giving the name of the added dimension in gf.

  • dim_units_in: A Unitful.Unitlike giving the units of the dimension as they appear in the NetCDF/MAT file.

  • dim_units_out: A String giving the kind of quantity being read; should be one of keys(UNITS).

Optional Arguments

  • transform: If provided, the dimension will be mapped according to transform before any other steps are taken. Default nothing.

  • force: If true, the range will be constructed even if the vector isn't linearly spaced by linear interpolation preserving the length. Default false.

source

Sargassum.add_temporal_dimension! Method
julia
add_temporal_dimension!(gf, infile, time_name_in, time_name_out, time_start, time_period; transform, force)

Add a new temporal dimension to gf::GriddedField with data read from a NetCDF or MAT file infile.

The new dimension appears last in the list of dimension names.

Arguments

  • gf: The GriddedField to be modified.

  • infile: The path to the NetCDF/MAT file.

  • time_name_in: A String giving the name of the time dimension to read in as it appears in the NetCDF/MAT file.

  • time_name_out: A Symbol giving the name of the added time dimension in gf.

  • time_start: A DateTime giving the reference time of the time dimension, e.g. if the units of the NetCDF/MAT are hours since 1990-01-01 then time_start == DateTime(1900, 1, 1).

  • time_period: A Unitful.FreeUnits giving the time step (units) of the time dimension. E.g. if the units of the NetCDF/MAT are hours since 1990-01-01 then time_period == u"hr".

Optional Arguments

  • transform: If provided, the dimension will be mapped according to transform before any other steps are taken. Default nothing.

  • force: If true, the range will be constructed even if the vector isn't linearly spaced by linear interpolation preserving the length. Default false.

source

Sargassum.ranges_increasing! Method
julia
ranges_increasing!(gf)

Moddify gf::GriddedField in place so that each dimension has variable which are increasing. Fields are automatically reversed if necessary.

source

Sargassum.sph2xy! Method
julia
sph2xy!(gf; lon_name, lat_name, x_name, y_name)

Moddify gf::GriddedField in place so that its longitudinal and latitudinal dimensions are converted to equirectangular coordinates.

Optional Arguments

  • lon_name: A Symbol giving the name of the longitudinal variable in gf.

  • lat_name: A Symbol giving the name of the latitudinal variable in gf.

  • x_name: A Symbol giving the name of the equirectangular x variable to be used in the modified gf.

  • y_name: A Symbol giving the name of the equirectangular y variable to be used in the modified gf.

source

Sargassum.LAND_ITP Constant
julia
const LAND_ITP

The interpolant for landmass location. This is a Ref, access or modify the actual interpolant with LAND_ITP.x.

source

Sargassum.NUTRIENTS_ITP Constant
julia
const NUTRIENTS_ITP

The interpolant for ocean nitrogen content. This is a Ref, access or modify the actual interpolant with NO3_ITP_ITP.x.

source

Sargassum.STOKES_ITP Constant
julia
const STOKES_ITP

The interpolant for Stokes drift velocity. This is a Ref, access or modify the actual interpolant with STOKES_ITP.x.

source

Sargassum.TEMPERATURE_ITP Constant
julia
const TEMPERATURE_ITP

The interpolant for ocean temperature. This is a Ref, access or modify the actual interpolant with TEMPERATURE_ITP.x.

source

Sargassum.WATER_ITP Constant
julia
const WATER_ITP

The interpolant for ocean currents. This is a Ref, access or modify the actual interpolant with WATER_ITP.x.

source

Sargassum.WAVES_ITP Constant
julia
const WAVES_ITP

The interpolant for wave height. This is a Ref, access or modify the actual interpolant with WAVES_ITP.x.

source

Sargassum.WIND_ITP Constant
julia
const WIND_ITP

The interpolant for wind speed. This is a Ref, access or modify the actual interpolant with WIND_ITP.x.

source

Sargassum.itps_load Method
julia
itps_load(; dir = _ITPS_SCRATCH.x)

Attempt to load the interpolants in directory dir.

This assumes that there exists in this directory the following files, each containing a variable as follows

File NameVariable Name
WATER_ITP.jld2WATER_ITP
WIND_ITP.jld2WIND_ITP
STOKES_ITP.jld2STOKES_ITP
WAVES_ITP.jld2WAVES_ITP
NUTRIENTS_ITP.jld2NUTRIENTS_ITP
TEMPERATURE_ITP.jld2TEMPERATURE_ITP
LAND_ITP.jld2LAND_ITP

source

Sargassum.itps_default_construct Method
julia
itps_default_construct(; download = false, verbose = true)

Construct all interpolants using the default data.

This overwrites any default interpolants already constructed.

Interpolants constructed: water, wind, stokes, waves, nutrients, temperature, land.

Optional Arguments

  • download: If true, download the data (roughly 1.2 GB of NetCDF files).

  • verbose: If true, print itp construction stats. Default true.

source

Sargassum.boundary Method
julia
boundary(itp)

Return a the corners of the spatio-temporal box the interpolant is defined in.

If itp is time-dependent, return (lon_min, lon_max, lat_min, lat_max, t_min, t_max).

If itp is time-independent, (lon_min, lon_max, lat_min, lat_max).

Assumes that the (x, y, t) variables are named (:x, :y, :t), respectively.

Examples

julia
julia> (lon_min, lon_max, lat_min, lat_max, t_min, t_max) = boundary(WATER_ITP)
julia
julia> (lon_min, lon_max, lat_min, lat_max) = boundary(LAND_ITP)

source

Sargassum.dim Method
julia
dim(itp, name)

Return the variable of itp whose value corresponds to the dimension indicated by name.

Use dims to see a list of possible values of name.

Example

julia
julia> x = dim(WATER_ITP, :x) # the x values defining interpolant knots

source

Sargassum.dims Method
julia
dims(itp)

Return the list of variable name/unit pairs of itp.

Example

julia
julia> dims(WATER_ITP)

source

Sargassum.field Method
julia
field(itp, name)

Return the sub-interpolant of itp whose value corresponds to the field indicated by name.

Use fields to see a list of possible values of name.

Example

julia
julia> v_x = field(WATER_ITP, :u) # the x component of the water velocity
julia> v_x(1, 2, 3) # evaluate it at `(x, y, t) = (1, 2, 3)`.

source

Sargassum.fields Method
julia
fields(itp)

Return the list of field name/unit pairs of itp.

Example

julia
julia> fields(WATER_ITP)

source

Sargassum.update_interpolant! Method
julia
update_interpolant!(itp, itp_new)

Update (replace) itp with itp_new, where itp_new should be an InterpolatedField and itp should be one of

  • WATER_ITP

  • WIND_ITP

  • STOKES_ITP

  • WAVES_ITP

  • NUTRIENTS_ITP

  • TEMPERATURE_ITP

  • LAND_ITP

source

RaftParameters

Rafts and Clumps

Sargassum.ClumpParameters Type
julia
struct ClumpParameters

A container for the high-level parameters of the BOM equations.

Fields

  • α []: The fraction of the wind field acting on the particle.

  • τ [d]: Measures the inertial response time of the medium to the particle

  • R []: A geometric parameter.

  • Ω [1/d]: The angular velocity of the Earth.

  • σ []: The Stokes drift parameter; this applies an additional fraction of the Stokes drift to the water velocity component of the particle.

Constructing from physical constants

ClumpParameters(; constants...)

Compute the parameters required for the eBOM equations from physical constants.

Constants

  • δ []: The bouyancy of the particle. Default: 1.25.

  • a [km]: The radius of the particle. Default: 1.0e-4.

  • ρ [kg/km^3]: The density of the water. Default: 1027.0e9.

  • ρa [kg/km^3]: The density of the air. Default: 1.2e9.

  • ν [km^2/d]: The viscosity of the water. Default: 8.64e-8.

  • νa [km^2/d]: The viscosity of the air. Default: 1.296e-6.

  • Ω [rad/d]: The angular velocity of the Earth. Default: .

  • σ []: The Stokes drift parameter. Default: 0.0.

Constructing directly

ClumpParameters(α, τ, R, Ω, σ)

Construct a ClumpParameters from the final constants directly. Use this to force specific values.

source

Sargassum.RaftParameters Type
julia
struct RaftParameters{S, C, G, L, I}

A container for the parameters defining a raft. Each clump and spring are identical.

Structure

RaftParameters acts as the parameter container for Raft!. The solution vector u is a 2 x N Matrix of the form [x1 x2 ... xN ; y1 y2 ... yN] giving the initial coordinates of each clump.

Fields

  • ics: An InitialConditions.

  • clumps: The ClumpParameters shared by each clump in the raft.

  • springs: A subtype of AbstractSpring.

  • connections: A subtype of AbstractConnections.

  • gd_model: A subtype of AbstractGrowthDeathModel.

  • land: A subtype of AbstractLand.

  • n_clumps_max: An Integer equal to the maximum allowed number of clumps. The number of clumps will not exceed this for any reason.

  • living: A BitVector such that living[i] == true if the clump with index i is alive.

  • n_clumps_tot: An Base.RefValue{Int64} whose reference is equal to the total number of clumps that have ever existed (i.e. it is at least the number of clumps that exist at any specific time.)

  • geometry: A Bool that toggles whether to apply the geometric correction factors γ_sphere and τ_sphere. Note that the simulation still uses the available interpolants, therefore if the interpolants have been created with geometric corrections included, but RaftParameters is created with geometry == false, the result will be a mixture of corrected and uncorrected terms.

  • dx_MR: dx of the Maxey-Riley equation. When provided, integration is done using FastRaft!.

  • dy_MR: dy of the Maxey-Riley equation. When provided, integration is done using FastRaft!.

Constructor

RaftParameters(; ics, clumps, springs, connections, gd_model, land, n_clumps_max, geometry = true, fast_raft = false)

The quantities living and n_clumps_tot are computed automatically under the assumption that the clumps initially provided are all alive.

Fast Raft

If fast_raft == true in the above constructor, then the equations will be integrated using FastRaft!. This is faster than Raft! at the expense of a more front-loaded computation since the interpolants must be computed. Using fast raft is advisable when the number of clumps is large. Default false.

One can also set fast_raft = (dx_MR, dy_MR) directly if the interpolants have been computed previously using dxdy_MR.

source

Sargassum.dxdy_MR Method
julia
dxdy_MR(tspan, clumps; geometry = true)

Compute (dx, dy) where dx and dy are interpolants evaluable at (x, y, t) equal to the right-hand-side of the Maxey-Riley equations (spring force excluded).

This is automatically applied when a fast raft is selected in Raft!.

Optional Arguments

  • geometry: Passed directly to γ_sphere](@ref) and τ_sphere.

source

Initial Conditions

Sargassum.InitialConditions Type
julia
struct InitialConditions

A container for the initial conditions for a raft.

Fields

  • tspan: A Tuple such that the integration is performed for tspan[1] ≤ t ≤ tspan[2] where t is either two two DateTimes or two Reals measured in UNITS["time"] since T_REF.

  • ics: A 2 x N Matrix of the form [x1 x2 ... xN ; y1 y2 ... yN] giving the initial coordinates of each clump.

Generic constructor

InitialConditions(;tspan, ics)

Constructing from positions

InitialConditions(tspan, xy0; to_xy)

Construct initial conditions suitable for use in RaftParameters.ics from 2 x N Matrix, xy0 which should be equirectangular coordinates.

InitialConditions(tspan, x_range, y_range; to_xy)

Generate clumps in a rectangular arrangement.

InitialConditions(tspan, x0, y0; to_xy)

Generate a single clump with coordinates (x0, y0).

Optional Arguments

to_xy: If true, the coordinates are converted from spherical to equirectangular coordinates. Default false.

Constructing from a SargassumDistribution

InitialConditions(tspan, dist, weeks, levels; seed)

Construct InitialConditions from a SargassumDistribution.

Arguments

  • tspan: A 2-Tuple such that the integration is performed for tspan[1] ≤ t ≤ tspan[2] where t is either two two DateTimes or two Reals measured in UNITS["time"] since T_REF.

  • dist: A SargassumDistribution.

  • weeks: A Vector{<:Integer} giving the weeks of the month to consider. Each entry should be between 1 and 4 and appear only once.

  • levels: The number of clump levels. Note that this is NOT equal to the number of clumps, see below.

Levels

Boxes with nonzero Sargassum content are divided into levels levels of size (maximum(D) - minimum(D))/levels where D = log10.(dist.sargassm[:,:,weeks]). Each box gets a number of clumps equal to its level index. For example, if levels = 2, then the smaller half of the boxes (by Sargassum content) get 1 clump each and the larger half get 2 clumps each.

Optional Arguments

  • seed: Random.seed!(seed) is called before the the initialization. Default 1234.

Constructing from a SargassumDistribution (streamlined)

InitialConditions(dist, week, n_weeks, levels; seed)

Construct InitialConditions from a SargassumDistribution using a more streamlined incantation than above.

Arguments

  • dist: A SargassumDistribution.

  • week: An integer between 1 and 4 giving the particular week of the month to initialize from.

  • n_weeks: The number of weeks to integrate for.

  • levels: As above.

Optional Arguments

  • seed: Random.seed!(seed) is called before the the initialization. Default 1234.

source

Springs

Sargassum.AbstractConnections Type
julia
abstract type AbstractConnections

A supertype for all connections between clumps.

Every subtype of AbstractConnections should have a field connections which is similar to a vector of vectors such that that connections[i] is a list of clump indices that are connected to clump i. connections should be initialized to vector of empty vectors of length n_clumps_max.

This should be updated in-place during the integration, i.e. it only shows the connections at the current time.

Every subtype of AbstractConnections should implement a form_connections(con::Connections, u) method which returns what con.connections should be updated with, assuming that u is the solution vector. The correction of indices due to living clumps is provided automatically later, so here it can be assumed that u contains only living clumps.

Any subtype of AbstractConnections can be evaluated at an OrdinaryDiffEq.integrator for callback purposes.

source

Sargassum.AbstractSpring Type
julia
abstract type AbstractSpring

A supertype for all spring parameters. Each clump, when conncted, is joined by the same kind of spring.

Every subtype of AbstractSpring should have a field k::Function representing the stiffness force and callable as k(x) as well as a field L::Real representing the spring's natural length.

All forces are computed using

julia
parameters.k(d)*(parameters.L/d - 1)*(xy1 - xy2)

where d = norm(xy1 - xy2).

source

Sargassum.BOMBSpring Type
julia
BOMBSpring{F}

A subtype of AbstractSpring representing a BOMB spring of the form A * (exp((x - 2*L)/0.2) + 1)^(-1).

Extra fields

  • A: The amplitude of the force.

Constructor

BOMBSpring(A::Real, L::Real)

source

Sargassum.ConnectionsFull Type
julia
struct ConnectionsFull

A connection type such that every clump is connected to every other clump.

Constructor

ConnectionsFull(n_clumps_max)

source

Sargassum.ConnectionsNearest Type
julia
struct ConnectionsNearest

A connection type such that every clump is connected to a number of its nearest neighbors.

Fields

  • neighbors: The number of nearest neighbors each clump should be connected to.

Constructors

ConnectionsNearest(n_clumps_max, neighbors)

source

Sargassum.ConnectionsNone Type
julia
struct ConnectionsNone

A connection type such that no clumps are connected.

Constructor

ConnectionsNone(n_clumps_max)

source

Sargassum.ConnectionsRadius Type
julia
struct ConnectionsRadius

A connection type such that every clump is connected to every clump within a given radius.

Fields

  • radius: A distance (assumed in UNITS["distance"]) such that each clump is connected to every clump whose distance is at most radius from it.

Constructor

ConnectionsRadius(n_clumps_max, radius)

source

Sargassum.HookeSpring Type
julia
HookeSpring{F}

A subtype of AbstractSpring representing a spring with a constant stiffness.

Constructor

HookeSpring(k::Real, L::Real)

source

Sargassum.ΔL Method
julia
ΔL(x_range, y_range; to_xy)

Compute a spring length from a rectangular arrangement of clumps provided by x_range and y_range. This is the distance between the centers of diagonally-adjacent gridpoints. These should be equirectangular coordinates.

Optional Arguments

to_xy: If true, the coordinates are converted from spherical to equirectangular coordinates. Default false.

source

Sargassum.ΔL Method
julia
ΔL(ics::InitialConditions; k::Integer)

Compute a spring length from a InitialConditions. This is the median among all pairwise equirectangular distances between points' k nearest neighbors. Default k = 5.

source

Sargassum.ΔL Method
julia
ΔL(dist::SargassumDistribution)

Compute a spring length from a SargassumDistribution. This is the equirectangular distance between the centers of diagonally-adjacent gridpoints.

source

Land

Sargassum.AbstractLand Type
julia
abstract type AbstractLand

A supertype for all land/shore types.

source

Sargassum.Land Type
julia
mutable struct Land{I}

A container for data handling death of clumps upon reaching the shore.

Fields

  • land_itp: A InterpolatedField such that land_itp.fields[:land](x, y) is equal to 1.0 if (x, y) is on land and 0.0 otherwise.

  • deaths: A Vector of indices of clumps that are to be killed.

  • verbose: A Bool such that verbose = true will log times and labels of clumps that hit land.

Constructors

Land(;land_itp::InterpolatedField = land_itp, verbose = false)

source

Sargassum.NoLand Type
julia
struct NoLand

An AbstractLand such that the land/shore is completely ignored.

source

Biology

Sargassum.AbstractGrowthDeathModel Type
julia
abstract type AbstractGrowthDeathModel

The abstract type for growth and death models.

source

Sargassum.BrooksModel Type
julia
mutable struct BrooksModel{B, D}

The growth/death model of Brooks et al. (2018).

Fields

  • S: A Vector{Float64} of length n_clumps_max representing an "amount" or "mass" for each clump.

  • S_gen. A Distributions.Sampleable{Univariate, ...} such that rand(S_gen) generates a sample of S. E.g. S_gen = Distributions.Dirac(0.0) will always initialize clumps with S = fill(0.0, n_clumps_max).

  • params: The BrooksModelParameters parameters of the model.

  • growths:A Vector of indices of clumps that are to be grown (if any).

  • deaths: A Vector of indices of clumps that are to be killed (if any).

  • verbose: A Bool such that verbose = true will log times and labels of clumps that grow and die.

Constructors

BrooksModel(n_clumps_max; S_gen = Dirac(0.0), params = BrooksModelParameters(), verbose = false)

Logic

At each time step, model.S is modified by model.params.dSdt * dt. The resulting S values are each compared to model.params.S_min and model.params.S_max and the associated ith clump either dies or spawns a child according to grow!(integrator, location = i).

source

Sargassum.BrooksModelParameters Type
julia
struct BrooksModelParameters{I, F}

A container for the parameters of the model of Brooks et al. (2018).

Parameters

  • temp [°C]: An InterpolatedField for the water temperature. Default TEMPERATURE_ITP.x.

  • no3 [mmol N/m^3]: An InterpolatedField for the Nitrogen content of the water. NUTRIENTS_ITP.x.

  • μ_max [1/d]: Sargassum maximum growth rate. Value: 0.1

  • m [1/d]: Sargassum mortality rate. Value: 0.05

  • k_N [mmol N/m^3]: Sargassum nutrient (N) uptake half saturation. Value: 0.012

  • T_min [°C]: Minimum temperature for Sargassum growth. Value: 10.0

  • T_max [°C]: Minimum temperature for Sargassum growth. Value: 40.0

  • clumps_limits: A Tuple of the form (n_clumps_min, n_clumps_max). These impose hard lower and upper limits on the total number of clumps that can exist at any specific time (the total number of clumps that can have ever existed - i.e. n_clumps_tot of RaftParameters - may be higher.) Default: (0, 10000).

  • S_min: A clump dies when S < S_min. Default 0.0.

  • S_max: A clump grows when S > S_max. Default 1.0.

  • dSdt: Compute the rate of change of the "amount" S according to the Brooks model.

dSdt

This function is of the form dSdt = growth_factors - death_factors.

  • growth_factors = μ_max * temperature_factor * nutrients_factor

  • death_factors = m

Constructors

BrooksModelParameters(; parameters...)

where each parameter is a kwarg with the default values given above.

source

Sargassum.ImmortalModel Type
julia
struct ImmortalModel

An AbstractGrowthDeathModel such that no growth or death occurs.

Constructors

ImmortalModel()

source

Sargassum.grow! Method
julia
grow!(integrator; location)

Add a clump to the RaftParameters, rp = integrator.p, with an index equal to rp.n_clumps_tot + 1 and also update rp.living appropriately.

Location

location can be a pre-defined flag, an integer, or a [x, y] vector. The default value is the flag "parent".

The possible flags are:

  • "parent": A parent clump is chosen randomly among clumps that already exist, and the new clump is placed a distance integrator.rp.springs.L away and at a uniformly random angle from it.

  • "com": The same as "parent", except the center location is at the center of mass of the raft.

If location is an Integer with value i, then the new clump will be grown with ith clump (by vector location) as its parent.

If location is a Vector{<:Real}, the new clump will be placed at those [x, y] coordinates.

source

Sargassum.kill! Method
julia
kill!(integrator, i)

Remove the clump with index i from integrator.

Can be applied as kill!(integrator, inds) in which case each i in inds will be killed in order.

For the RaftParameters, rp = integrator.p, update rp.living appropriately.

source

Simulation

Sargassum.FastRaft! Method
julia
FastRaft!(du, u, p, t)

When integrated, produces a result (nearly) identical to Raft!, but is generally faster at the expense of a more front-loaded computation due to the requirement of additional interpolants.

source

Sargassum.Leeway! Method
julia
Leeway!(du, u, p::RaftParameters, t)

Compute the right-hand-side of the differential equation controlling the motion of raft particles whose velocities are equal to u = v_water + α v_wind.

The parameters p are given by RaftParameters, but only p.α is used.

source

Sargassum.Raft! Method
julia
Raft!(du, u, p::RaftParameters, t)

Compute the right-hand-side of the differential equation controlling the motion of a raft with parameters given by RaftParameters.

The solution vector u is a 2 x N Matrix of the form [x1 x2 ... xN ; y1 y2 ... yN] giving the coordinates of each clump.

For integrating using a leeway velocity, Leeway! should be used.

source

Sargassum.simulate Method
julia
simulate(rp; leeway, alg, abstol, reltol, showprogress, dt, return_raw)

Simulate a Sargassum raft with RaftParameters rp and return a RaftTrajectory.

This function modifies the fields of rp significantly.

Arguments

Optional Arguments

  • leeway: Use Leeway! to integrate particles with no springs or inertia. Default false.

  • alg: The integration algorithm to use, default Tsit5().

  • abstol: The absolute tolerance of integration; default nothing.

  • reltol: The relative tolerance of integration; default nothing.

  • showprogress: If true, print a status indicator of the progress of the integration. Default false.

  • dt: The solution trajectories are uniformized to be spaced in time by increments of dt. Note that the units of this quantity are implicity UNITS["time"]. Default 0.1.

  • return_raw: If true, return the result of OrdinaryDiffEq.solve, rather than a RaftTrajectory. Use this if you would like to manipulate the solution directly. Default false.

source

Sargassum.RaftTrajectory Type
julia
struct RaftTrajectory

A container for the data of a every clump's trajectory in a raft, as well as its center of mass.

Fields

  • trajectories: A Dict mapping clump indices to their corresponding Trajectory.

  • t: A vector of all time possible slices across the clump trajectories.

  • n_clumps: A vector such that n_clumps[i] is the number of clumps that are alive at time t[i].

  • com: A Trajectory corresponding to the center of mass of the raft.

Constructor

RaftTrajectory(; trajectories, n_clumps, com)

The field t is set to com.t.

Plotting

This object can be viz and viz!.

source

Sargassum.Trajectory Type
julia
struct Trajectory

A container for the data of a single clump's trajectory.

Fields

  • xy: A Matrix of size N x 2 such that xy[i,:] gives the [x, y] or [lon, lat] coordinates at the clump at time t[i].

  • t: A Vector of length N giving the time values of the trajectory.

Constructor

Trajectory(xy, t)

Construct a Trajectory directly from xy and t as defined above.

Plotting

This object can be viz and viz!.

source

Sargassum.bins Method
julia
bins(raft_trajectory, dist; return_xy_bins)

Equivalent to bins(raft_trajectory, x_bins, y_bins where x_bins and y_bins are computed automatically from the SargassumDistribution, dist.lon and dist.lat.

This assumes that dist.lon and dist.lat give the central locations of the dist bins.

Optional Arguments

  • return_xy_bins: A Bool such that, if true, the tuple (x_bins, y_bins, bins) is returned instead of just bins. Default false.

source

Sargassum.bins Method
julia
bins(raft_trajectory, x_bins, y_bins)

Return a matrix mat such that mat[i, j] is the number of points in raft_trajectory that, at any time, were inside the rectangle lon ∈ (x_bins[i], x_bins[i + 1]), lat ∈ (y_bins[i], y_bins[i + 1]).

Both x_bins and y_bins should be StepRangeLen, i.e. of the form range(start, stop, length = L). Then, mat has dimensions length(x_bins) - 1 x length(y_bins) - 1.

No coversion from or to spherical coordinates is done on x_bins and y_bins.

source

Sargassum.time_slice Method
julia
time_slice(traj, tspan)

Return a new Trajectory consisting of points and times of traj that are between first(tspan) and last(tspan). The result Trajectory may be empty.

Can also be applied to a RaftTrajectory in which case time_slice is applied to each member Trajectory.

source

I/O

Sargassum.rtr2mat Method
julia
rtr2mat(rtr, outfile; force)

Write the RaftTrajectory in rtr to outfile which must be a .mat file.

This writes the raw, unbinned trajectory data.

Optional Arguments

  • force: If true, delete outfile if it already exists. Default false.

source

Sargassum.rtr2nc Method
julia
rtr2nc(rtr, outfile, dist; force)

Write the RaftTrajectory in rtr to outfile which must be a .nc file.

The data are binned by passing dist to bins, i.e. the bins are chosen to be the same as the bins of dist.

It is required that rtr.t is linearly spaced.

Optional Arguments

  • force: If true, delete outfile if it already exists. Default false.

source

Sargassum.rtr2nc Method
julia
rtr2nc(rtr, outfile, lon_bins, lat_bins; force)

Write the RaftTrajectory in rtr to outfile which must be a .nc file.

The data are binned by passing lon_bins and lat_bins to bins.

It is required that rtr.t is linearly spaced.

Optional Arguments

  • force: If true, delete outfile if it already exists. Default false.

source

Examples

Sargassum.QuickRaftParameters Method
julia
QuickRaftParameters(ics; kwargs...)

Generate a RaftParameters object to integrate from ics::InitialConditions.

Optional Arguments

  • use_biology: If true, include biological effects in the simulation. Default false.

  • n_connections: The number of inter-clump nearest neighbor connectionsto form. Default 2.

  • delta: The bouyancy of the particle. Default: 1.14.

  • tau: The inertial response time, measured in days. Default 0.0103.

  • sigma: The Stokes drift coefficient. Default 1.20.

  • A_spring: The magnitude of the eBOMB spring stiffness. Default 15.1.

  • lambda_spring: A factor multiplying the springs' natural lengths. Default 2.97.

  • mu_max: Sargassum maximum growth rate, measured in inverse days. Default 0.00542.

  • m: Sargassum mortality rate, measured in inverse days. Default 0.00403.

  • k_N: Sargassum nutrient (N) uptake half saturation, measured in mmol N/m^3. Default 0.000129.

  • S_min: A clump dies when it's "amount" drops below this value. Default -0.00481.

  • S_max: A clump dies when it's "amount" grows above this value. Default 0.001.

  • geometry: A Bool that toggles whether to apply the geometric correction factors γ_sphere and τ_sphere. Note that the simulation still uses the available interpolants, therefore if the interpolants have been created with geometric corrections included, but RaftParameters is created with geometry == false, the result will be a mixture of corrected and uncorrected terms. Default true.

  • verbose: Whether to print out clump growths/deaths at each step. Default false.

source

Sargassum.QuickRaftParameters Method
julia
QuickRaftParameters()

Return a simple RaftParameters with fixed parameters suitable for testing purposes.

source