Simulation API
These are the full docstrings for the Simulation subsection of Sargassum.jl.
Index
Sargassum.LAND_ITPSargassum.NUTRIENTS_ITPSargassum.STOKES_ITPSargassum.TEMPERATURE_ITPSargassum.WATER_ITPSargassum.WAVES_ITPSargassum.WIND_ITPSargassum.AbstractConnectionsSargassum.AbstractGrowthDeathModelSargassum.AbstractLandSargassum.AbstractSpringSargassum.BOMBSpringSargassum.BrooksModelSargassum.BrooksModelParametersSargassum.ClumpParametersSargassum.ConnectionsFullSargassum.ConnectionsNearestSargassum.ConnectionsNoneSargassum.ConnectionsRadiusSargassum.GriddedFieldSargassum.HookeSpringSargassum.ImmortalModelSargassum.InitialConditionsSargassum.InterpolatedFieldSargassum.LandSargassum.NoLandSargassum.RaftParametersSargassum.RaftTrajectorySargassum.TrajectorySargassum.FastRaft!Sargassum.Leeway!Sargassum.QuickRaftParametersSargassum.QuickRaftParametersSargassum.Raft!Sargassum.add_derivatives!Sargassum.add_field!Sargassum.add_spatial_dimension!Sargassum.add_temporal_dimension!Sargassum.binsSargassum.binsSargassum.boundarySargassum.dimSargassum.dimsSargassum.dxdy_MRSargassum.fieldSargassum.fieldsSargassum.grow!Sargassum.itps_default_constructSargassum.itps_loadSargassum.kill!Sargassum.ranges_increasing!Sargassum.rtr2matSargassum.rtr2ncSargassum.rtr2ncSargassum.simulateSargassum.sph2xy!Sargassum.time_sliceSargassum.update_interpolant!Sargassum.ΔLSargassum.ΔLSargassum.ΔL
Interpolants
Sargassum.GriddedField Type
struct GriddedField{N, T, U}A container for gridded data, possibly time-dependent.
Fields
dims_names: AVectorofTuple{Symbol, Unitful.Unitlike}s such thatdims_names[i][1]is theith dimension offieldsanddims_names[i][2]gives the units of theith dimension.dims: ADictmapping variable names to ranges they take.fields_names: AVectorofTuple{Symbol, Unitful.Unitlike}s such thatfields_names[i][1]is theith field andfields_names[i][2]gives its units.fields: ADictmapping 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.
Sargassum.InterpolatedField Type
struct InterpolatedField{N, T, U, I}A container for interpolants of gridded data, possibly time-dependent.
Fields
dims_names: AVectorofTuple{Symbol, Unitful.Unitlike}s such thatdims_names[i][1]is theith dimension offieldsanddims_names[i][2]gives the units of theith dimension.dims: ADictmapping variable names to ranges they take.fields_names: AVectorofTuple{Symbol, Unitful.Unitlike}s such thatfields_names[i][1]is theith field andfields_names[i][2]gives its units.fields: ADictmapping 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, anyInterpolations.InterpolationTypecan be provided. Default"cubic".extrapolate_value: A constant extrapolation is performed with this value. Default"0.0".
Plotting
Sargassum.add_derivatives! Method
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: AnInterpolatedField. The field should havex,yandtvariables along withxandycomponents 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, anyInterpolations.InterpolationTypecan be provided. Default"cubic".extrapolate_value: A constant extrapolation is performed with this value. Default"0.0".xyt_names: ATuplewith three symbols, corresponding to thex,yandtvariables in that order. Default(:x, :y, :t).vxvy_names: ATuplewith two symbols, corresponding to thexandycomponents of the vector field, in that order. Default(:u, :v).Dx_Dy_vort_names: ATuplewith three symbols, corresponding to thexcomponent of the material derivative, theycomponent of the material derivative and the vorticity, in that order. Default(:DDt_x, :DDt_y, :vorticity).geometry: Iftrue, include factors that take into account the spherical geometry of the Earth. Defaulttrue.
Sargassum.add_field! Method
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: TheGriddedFieldto be modified.infile: The path to the NetCDF/MAT file.field_name_in: AStringgiving the name of the field to read in as it appears in the NetCDF/MAT file.field_name_out: ASymbolgiving the name of the added field ingf.field_units_in: AUnitful.Unitlikegiving the units of the field as they appear in the NetCDF/MAT file.field_units_out: AStringgiving the kind of quantity being read; should be one ofkeys(UNITS).
Optional Arguments
take_axes: If provided, only the selected elements will be taken viafield[take_axes...]. For example, if thefieldis four dimensional, passingtake_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. Defaultnothing.permutation: If provided, the field will be permuted according topermutation. Applied AFTERtake_axes. Defaultnothing.scale_factor_name: The name of the scale factor, only for NetCDF files. If no scale factor is found, it is taken to be1. 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 be0. Default"add_offset".missings_name: A vector of names of missing/fill/extra values, only for NetCDF files. Each such value will be replaced bymissings_replacementif found. Default["_FillValue", "missing_value"].missings_replacement:missings_namereplaces the missing/fill/extra values with this. Default0.0.
Sargassum.add_spatial_dimension! Method
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: TheGriddedFieldto be modified.infile: The path to the NetCDF/MAT file.dim_name_in: AStringgiving the name of the dimension to read in as it appears in the NetCDF/MAT file.dim_name_out: ASymbolgiving the name of the added dimension ingf.dim_units_in: AUnitful.Unitlikegiving the units of the dimension as they appear in the NetCDF/MAT file.dim_units_out: AStringgiving the kind of quantity being read; should be one ofkeys(UNITS).
Optional Arguments
transform: If provided, the dimension will be mapped according totransformbefore any other steps are taken. Defaultnothing.force: Iftrue, the range will be constructed even if the vector isn't linearly spaced by linear interpolation preserving the length. Defaultfalse.
Sargassum.add_temporal_dimension! Method
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: TheGriddedFieldto be modified.infile: The path to the NetCDF/MAT file.time_name_in: AStringgiving the name of the time dimension to read in as it appears in the NetCDF/MAT file.time_name_out: ASymbolgiving the name of the added time dimension ingf.time_start: ADateTimegiving the reference time of the time dimension, e.g. if the units of the NetCDF/MAT arehours since 1990-01-01thentime_start == DateTime(1900, 1, 1).time_period: AUnitful.FreeUnitsgiving the time step (units) of the time dimension. E.g. if the units of the NetCDF/MAT arehours since 1990-01-01thentime_period == u"hr".
Optional Arguments
transform: If provided, the dimension will be mapped according totransformbefore any other steps are taken. Defaultnothing.force: Iftrue, the range will be constructed even if the vector isn't linearly spaced by linear interpolation preserving the length. Defaultfalse.
Sargassum.ranges_increasing! Method
ranges_increasing!(gf)Moddify gf::GriddedField in place so that each dimension has variable which are increasing. Fields are automatically reversed if necessary.
Sargassum.sph2xy! Method
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: ASymbolgiving the name of the longitudinal variable ingf.lat_name: ASymbolgiving the name of the latitudinal variable ingf.x_name: ASymbolgiving the name of the equirectangularxvariable to be used in the modifiedgf.y_name: ASymbolgiving the name of the equirectangularyvariable to be used in the modifiedgf.
Sargassum.LAND_ITP Constant
const LAND_ITPThe interpolant for landmass location. This is a Ref, access or modify the actual interpolant with LAND_ITP.x.
Sargassum.NUTRIENTS_ITP Constant
const NUTRIENTS_ITPThe interpolant for ocean nitrogen content. This is a Ref, access or modify the actual interpolant with NO3_ITP_ITP.x.
Sargassum.STOKES_ITP Constant
const STOKES_ITPThe interpolant for Stokes drift velocity. This is a Ref, access or modify the actual interpolant with STOKES_ITP.x.
Sargassum.TEMPERATURE_ITP Constant
const TEMPERATURE_ITPThe interpolant for ocean temperature. This is a Ref, access or modify the actual interpolant with TEMPERATURE_ITP.x.
Sargassum.WATER_ITP Constant
const WATER_ITPThe interpolant for ocean currents. This is a Ref, access or modify the actual interpolant with WATER_ITP.x.
Sargassum.WAVES_ITP Constant
const WAVES_ITPThe interpolant for wave height. This is a Ref, access or modify the actual interpolant with WAVES_ITP.x.
Sargassum.WIND_ITP Constant
const WIND_ITPThe interpolant for wind speed. This is a Ref, access or modify the actual interpolant with WIND_ITP.x.
Sargassum.itps_load Method
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 Name | Variable Name |
|---|---|
| WATER_ITP.jld2 | WATER_ITP |
| WIND_ITP.jld2 | WIND_ITP |
| STOKES_ITP.jld2 | STOKES_ITP |
| WAVES_ITP.jld2 | WAVES_ITP |
| NUTRIENTS_ITP.jld2 | NUTRIENTS_ITP |
| TEMPERATURE_ITP.jld2 | TEMPERATURE_ITP |
| LAND_ITP.jld2 | LAND_ITP |
Sargassum.itps_default_construct Method
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: Iftrue, download the data (roughly 1.2 GB of NetCDF files).verbose: Iftrue, print itp construction stats. Defaulttrue.
Sargassum.boundary Method
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> (lon_min, lon_max, lat_min, lat_max, t_min, t_max) = boundary(WATER_ITP)julia> (lon_min, lon_max, lat_min, lat_max) = boundary(LAND_ITP)Sargassum.dim Method
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> x = dim(WATER_ITP, :x) # the x values defining interpolant knotsSargassum.dims Method
dims(itp)Return the list of variable name/unit pairs of itp.
Example
julia> dims(WATER_ITP)Sargassum.field Method
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> 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)`.Sargassum.fields Method
fields(itp)Return the list of field name/unit pairs of itp.
Example
julia> fields(WATER_ITP)Sargassum.update_interpolant! Method
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_ITPWIND_ITPSTOKES_ITPWAVES_ITPNUTRIENTS_ITPTEMPERATURE_ITPLAND_ITP
RaftParameters
Rafts and Clumps
Sargassum.ClumpParameters Type
struct ClumpParametersA 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 particleR[]: 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:2π.σ[]: 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.
Sargassum.RaftParameters Type
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: AnInitialConditions.clumps: TheClumpParametersshared by each clump in the raft.springs: A subtype ofAbstractSpring.connections: A subtype ofAbstractConnections.gd_model: A subtype ofAbstractGrowthDeathModel.land: A subtype ofAbstractLand.n_clumps_max: AnIntegerequal to the maximum allowed number of clumps. The number of clumps will not exceed this for any reason.living: ABitVectorsuch thatliving[i] == trueif the clump with indexiis alive.n_clumps_tot: AnBase.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: ABoolthat toggles whether to apply the geometric correction factorsγ_sphereandτ_sphere. Note that the simulation still uses the available interpolants, therefore if the interpolants have been created with geometric corrections included, butRaftParametersis created withgeometry == false, the result will be a mixture of corrected and uncorrected terms.dx_MR:dxof the Maxey-Riley equation. When provided, integration is done usingFastRaft!.dy_MR:dyof the Maxey-Riley equation. When provided, integration is done usingFastRaft!.
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.
Sargassum.dxdy_MR Method
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.
Initial Conditions
Sargassum.InitialConditions Type
struct InitialConditionsA container for the initial conditions for a raft.
Fields
tspan: ATuplesuch that the integration is performed fortspan[1] ≤ t ≤ tspan[2]wheretis either two twoDateTimes or twoReals measured inUNITS["time"]sinceT_REF.ics: A2 x NMatrixof 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-Tuplesuch that the integration is performed fortspan[1] ≤ t ≤ tspan[2]wheretis either two twoDateTimes or twoReals measured inUNITS["time"]sinceT_REF.dist: ASargassumDistribution.weeks: AVector{<: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: ASargassumDistribution.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.
Springs
Sargassum.AbstractConnections Type
abstract type AbstractConnectionsA 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.
Sargassum.AbstractSpring Type
abstract type AbstractSpringA 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
parameters.k(d)*(parameters.L/d - 1)*(xy1 - xy2)where d = norm(xy1 - xy2).
Sargassum.BOMBSpring Type
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)Sargassum.ConnectionsFull Type
struct ConnectionsFullA connection type such that every clump is connected to every other clump.
Constructor
ConnectionsFull(n_clumps_max)Sargassum.ConnectionsNearest Type
struct ConnectionsNearestA 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)Sargassum.ConnectionsNone Type
struct ConnectionsNoneA connection type such that no clumps are connected.
Constructor
ConnectionsNone(n_clumps_max)Sargassum.ConnectionsRadius Type
struct ConnectionsRadiusA connection type such that every clump is connected to every clump within a given radius.
Fields
radius: A distance (assumed inUNITS["distance"]) such that each clump is connected to every clump whose distance is at mostradiusfrom it.
Constructor
ConnectionsRadius(n_clumps_max, radius)Sargassum.HookeSpring Type
HookeSpring{F}A subtype of AbstractSpring representing a spring with a constant stiffness.
Constructor
HookeSpring(k::Real, L::Real)Sargassum.ΔL Method
Δ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.
Sargassum.ΔL Method
Δ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.
Sargassum.ΔL Method
ΔL(dist::SargassumDistribution)Compute a spring length from a SargassumDistribution. This is the equirectangular distance between the centers of diagonally-adjacent gridpoints.
Land
Sargassum.Land Type
mutable struct Land{I}A container for data handling death of clumps upon reaching the shore.
Fields
land_itp: AInterpolatedFieldsuch thatland_itp.fields[:land](x, y)is equal to1.0if(x, y)is on land and0.0otherwise.deaths: AVectorof indices of clumps that are to be killed.verbose: ABoolsuch thatverbose = truewill log times and labels of clumps that hit land.
Constructors
Land(;land_itp::InterpolatedField = land_itp, verbose = false)Sargassum.NoLand Type
struct NoLandAn AbstractLand such that the land/shore is completely ignored.
Biology
Sargassum.AbstractGrowthDeathModel Type
abstract type AbstractGrowthDeathModelThe abstract type for growth and death models.
Sargassum.BrooksModel Type
mutable struct BrooksModel{B, D}The growth/death model of Brooks et al. (2018).
Fields
S: AVector{Float64}of lengthn_clumps_maxrepresenting an "amount" or "mass" for each clump.S_gen. ADistributions.Sampleable{Univariate, ...}such thatrand(S_gen)generates a sample ofS. E.g.S_gen = Distributions.Dirac(0.0)will always initialize clumps withS = fill(0.0, n_clumps_max).params: TheBrooksModelParametersparameters of the model.growths:AVectorof indices of clumps that are to be grown (if any).deaths: AVectorof indices of clumps that are to be killed (if any).verbose: ABoolsuch thatverbose = truewill 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).
Sargassum.BrooksModelParameters Type
struct BrooksModelParameters{I, F}A container for the parameters of the model of Brooks et al. (2018).
Parameters
temp[°C]: AnInterpolatedFieldfor the water temperature. DefaultTEMPERATURE_ITP.x.no3[mmol N/m^3]: AnInterpolatedFieldfor the Nitrogen content of the water.NUTRIENTS_ITP.x.μ_max[1/d]: Sargassum maximum growth rate. Value:0.1m[1/d]: Sargassum mortality rate. Value:0.05k_N[mmol N/m^3]: Sargassum nutrient (N) uptake half saturation. Value:0.012T_min[°C]: Minimum temperature for Sargassum growth. Value:10.0T_max[°C]: Minimum temperature for Sargassum growth. Value:40.0clumps_limits: ATupleof 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_totofRaftParameters- may be higher.) Default:(0, 10000).S_min: A clump dies whenS < S_min. Default0.0.S_max: A clump grows whenS > S_max. Default1.0.dSdt: Compute the rate of change of the "amount"Saccording to the Brooks model.
dSdt
This function is of the form dSdt = growth_factors - death_factors.
growth_factors = μ_max * temperature_factor * nutrients_factordeath_factors = m
Constructors
BrooksModelParameters(; parameters...)where each parameter is a kwarg with the default values given above.
Sargassum.ImmortalModel Type
struct ImmortalModelAn AbstractGrowthDeathModel such that no growth or death occurs.
Constructors
ImmortalModel()Sargassum.grow! Method
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 distanceintegrator.rp.springs.Laway 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.
Sargassum.kill! Method
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.
Simulation
Sargassum.FastRaft! Method
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.
Sargassum.Leeway! Method
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.
Sargassum.Raft! Method
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.
Sargassum.simulate Method
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
rp: ARaftParametersdefining the raft.
Optional Arguments
leeway: UseLeeway!to integrate particles with no springs or inertia. Defaultfalse.alg: The integration algorithm to use, defaultTsit5().abstol: The absolute tolerance of integration; defaultnothing.reltol: The relative tolerance of integration; defaultnothing.showprogress: Iftrue, print a status indicator of the progress of the integration. Defaultfalse.dt: The solution trajectories are uniformized to be spaced in time by increments ofdt. Note that the units of this quantity are implicityUNITS["time"]. Default0.1.return_raw: If true, return the result ofOrdinaryDiffEq.solve, rather than aRaftTrajectory. Use this if you would like to manipulate the solution directly. Defaultfalse.
Sargassum.RaftTrajectory Type
struct RaftTrajectoryA container for the data of a every clump's trajectory in a raft, as well as its center of mass.
Fields
trajectories: ADictmapping clump indices to their correspondingTrajectory.t: A vector of all time possible slices across the clump trajectories.n_clumps: A vector such thatn_clumps[i]is the number of clumps that are alive at timet[i].com: ATrajectorycorresponding to the center of mass of the raft.
Constructor
RaftTrajectory(; trajectories, n_clumps, com)The field t is set to com.t.
Plotting
Sargassum.Trajectory Type
struct TrajectoryA container for the data of a single clump's trajectory.
Fields
xy: AMatrixof sizeN x 2such thatxy[i,:]gives the[x, y]or[lon, lat]coordinates at the clump at timet[i].t: AVectorof lengthNgiving the time values of the trajectory.
Constructor
Trajectory(xy, t)Construct a Trajectory directly from xy and t as defined above.
Plotting
Sargassum.bins Method
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: ABoolsuch that, iftrue, the tuple(x_bins, y_bins, bins)is returned instead of justbins. Defaultfalse.
Sargassum.bins Method
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.
Sargassum.time_slice Method
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.
I/O
Sargassum.rtr2mat Method
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: Iftrue, deleteoutfileif it already exists. Defaultfalse.
Sargassum.rtr2nc Method
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: Iftrue, deleteoutfileif it already exists. Defaultfalse.
Sargassum.rtr2nc Method
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: Iftrue, deleteoutfileif it already exists. Defaultfalse.
Examples
Sargassum.QuickRaftParameters Method
QuickRaftParameters(ics; kwargs...)Generate a RaftParameters object to integrate from ics::InitialConditions.
Optional Arguments
use_biology: Iftrue, include biological effects in the simulation. Defaultfalse.n_connections: The number of inter-clump nearest neighbor connectionsto form. Default2.delta: The bouyancy of the particle. Default:1.14.tau: The inertial response time, measured in days. Default0.0103.sigma: The Stokes drift coefficient. Default1.20.A_spring: The magnitude of the eBOMB spring stiffness. Default15.1.lambda_spring: A factor multiplying the springs' natural lengths. Default2.97.mu_max: Sargassum maximum growth rate, measured in inverse days. Default0.00542.m: Sargassum mortality rate, measured in inverse days. Default0.00403.k_N: Sargassum nutrient (N) uptake half saturation, measured in mmol N/m^3. Default0.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. Default0.001.geometry: ABoolthat toggles whether to apply the geometric correction factorsγ_sphereandτ_sphere. Note that the simulation still uses the available interpolants, therefore if the interpolants have been created with geometric corrections included, butRaftParametersis created withgeometry == false, the result will be a mixture of corrected and uncorrected terms. Defaulttrue.verbose: Whether to print out clump growths/deaths at each step. Defaultfalse.
Sargassum.QuickRaftParameters Method
QuickRaftParameters()Return a simple RaftParameters with fixed parameters suitable for testing purposes.