Ideal age
The tracer equation for the ideal age is
\[\left(\partial_t + \mathbf{T}\right) \boldsymbol{a} = 1 - \frac{\boldsymbol{a}}{τ} \, (\boldsymbol{z} \le z_0),\]
where the sink term on the right clamps the age to $0$ at the surface (where $\boldsymbol{z} \le z_0$). The smaller the timescale $\tau$, the quicker $\boldsymbol{a}$ is restored to $0$ at the surface.
AIBECS can interpret tracer equations as long as you arrange them under the generic form:
\[\big(\partial_t + \mathbf{T}(\boldsymbol{p}) \big) \boldsymbol{x} = \boldsymbol{G}(\boldsymbol{x}, \boldsymbol{p}),\]
where $\mathbf{T}(\boldsymbol{p})$ is the transport, $\boldsymbol{G}(\boldsymbol{x}, \boldsymbol{p})$ is the net local sources and sinks, and $\boldsymbol{p}$ is the vector of model parameters. We will then use the AIBECS to simulate the ideal age by finding the steady-state of the system, i.e., the solution of
\[\partial_t \boldsymbol{x} = \boldsymbol{F}(\boldsymbol{x}, \boldsymbol{p}) = \boldsymbol{G}(\boldsymbol{x}, \boldsymbol{p}) - \mathbf{T}(\boldsymbol{p}) \, \boldsymbol{x} = 0.\]
In this tutorial, we will simulate the ideal age by
- defining functions for
T(p)andG(x,p), - defining the parameters
p, - generating the state function
F(x,p)and solving the associated steady-state problem, - and finally making a plot of our simulated ideal age.
We start by telling Julia that we want to use the AIBECS package and the OCIM2 circulation (the Ocean Circulation Inverse Model[1]).
using AIBECS
grd, TOCIM2 = OCIM2.load()(,
⠻⣦⡙⢦⡙⢦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠳⣌⠻⣦⡙⢦⡙⢦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠳⣌⠳⣌⠻⣦⡙⢦⡙⢦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠈⠳⣌⠳⣌⠻⣦⡙⢦⡙⢦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠈⠳⣌⠳⣌⠻⣦⡙⢦⡙⢦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠈⠳⣌⠳⣌⠻⣦⡙⣦⡙⣆⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠈⠳⣌⠳⣬⡻⣮⡻⣮⡳⣄⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠳⢬⡻⣮⡻⣮⡳⣌⠳⣄⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢮⡙⢮⡻⣮⡳⣌⠳⣄⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡙⢮⡻⣮⡳⣌⠳⣄⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡙⢮⡻⣮⡳⣌⠳⣄⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡙⢮⡻⣮⡳⣌⠳⣄⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡙⢮⡻⣮⡳⣌⠳⣄⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡙⢮⡻⣮⡳⣌⠳⣄⠀⠀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡙⢮⡻⣮⡳⣌⢳⡀⠀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡙⢮⡻⣮⡳⣝⢦⡀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⠲⣝⢮⡻⣮⡳⣝⢦⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠳⣝⢮⡻⣮⣯⡳⡄⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠓⢯⡻⣿⣿⣿⡆
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠉⠻⠿⣿⣿)If it's your first time, Julia will ask you to download the OCIM2, in which case you should accept (i.e., type y and "return"). Once downloaded, AIBECS will remember where it downloaded the file and it will only load it from your laptop.
grd is an OceanGrid object containing information about the 3D grid of the OCIM2 circulation and TOCIM2 is the transport matrix representing advection and diffusion.
The local sources and sinks for the age take the form
function G(x,p)
@unpack τ, z₀ = p
return @. 1 - x / τ * (z ≤ z₀)
endG (generic function with 1 method)as per the tracer equation. The @unpack line unpacks the parameters τ and z₀. The return line returns the net sources and sinks. (The @. "macro" tells Julia that the operations apply to every element.)
We can define the vector z of depths with depthvec.
z = depthvec(grd)200160-element Vector{Float64}:
18.0675569520817
18.0675569520817
18.0675569520817
18.0675569520817
18.0675569520817
18.0675569520817
18.0675569520817
18.0675569520817
18.0675569520817
18.0675569520817
⋮
5433.2531421838175
5433.2531421838175
5433.2531421838175
5433.2531421838175
5433.2531421838175
5433.2531421838175
5433.2531421838175
5433.2531421838175
5433.2531421838175Now we must construct a type for p the parameters. This type must contain our parameters τ and z₀.
struct IdealAgeParameters{U} <: AbstractParameters{U}
τ::U
z₀::U
endThe type is now ready for us to generate an instance of the parameter p. Let's use τ = 1.0 (s) and z₀ the minimum depth of the model.
p = IdealAgeParameters(1.0, 30.0) Row │ Symbol Value
│ Symbol Float64
─────┼─────────────────
1 │ τ 1.0
2 │ z₀ 30.0We now use the AIBECS to generate the state function $\boldsymbol{F}$ (and its Jacobian) via
F = AIBECSFunction(TOCIM2, G)(::SciMLBase.ODEFunction{false, AIBECS.var"#f#71"{Tuple{AIBECS.var"#49#50"{SparseMatrixCSC{Float64, Int64}}}, Vector{Int64}, AIBECS.var"#G#69"{Tuple{typeof(Main.G)}, AIBECS.var"#tracers#67"{Int64, Int64}}, AIBECS.var"#tracer#68"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#76"{AIBECS.var"#T#73"{Tuple{AIBECS.var"#49#50"{SparseMatrixCSC{Float64, Int64}}}, Int64, Vector{Int64}}, AIBECS.var"#∇ₓG#72"{Tuple{typeof(Main.G)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}) (generic function with 7 methods)Now that F and p are defined, we are going to solve for the steady-state. But first, we must create a SteadyStateProblem object that contains F, p, and an initial guess x_init for the age. (SteadyStateProblem comes from DiffEqBase.)
Let's make a vector of 0's for our initial guess.
nb = sum(iswet(grd)) # number of wet boxes
x_init = zeros(nb) # Start with age = 0 everywhere200160-element Vector{Float64}:
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
⋮
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0Now we can create our SteadyStateProblem instance
prob = SteadyStateProblem(F, x_init, p)SteadyStateProblem with uType Vector{Float64}. In-place: false
u0: 200160-element Vector{Float64}:
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
⋮
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0And finally, we can solve this problem, using the AIBECS CTKAlg() algorithm,
age = solve(prob, CTKAlg())u: 200160-element Vector{Float64}:
804.9793529617989
2704.2854677669366
559.4624759776623
277.68101763734074
571.2245067452855
406.14915296475453
416.88383377410975
645.3785433171208
411.6472087934807
242.71573068291764
⋮
1.2204098653557419e10
1.182084183210928e10
1.378396033391721e10
1.3083180539737455e10
1.1999572353055992e10
1.1227617212689003e10
1.3402055608211182e10
1.33226551540552e10
1.2121012215573372e10This should take a few seconds.
To conclude this tutorial, let's have a look at the age using AIBECS' plotting recipes and Plots.jl.
using PlotsWe first convert the age in years (because the default SI unit we used, i.e., seconds, is a bit small relative to global ocean timescales).
age_in_yrs = age * u"s" .|> u"yr"200160-element Vector{Quantity{Float64, 𝐓, Unitful.FreeUnits{(yr,), 𝐓, nothing}}}:
2.5508256425133687e-5 yr
8.569363537680104e-5 yr
1.7728296067434224e-5 yr
8.799180471180975e-6 yr
1.8101012331269976e-5 yr
1.2870090024740618e-5 yr
1.3210251532883036e-5 yr
2.045081195392301e-5 yr
1.3044312900647725e-5 yr
7.691197387726494e-6 yr
⋮
386.72454982499994 yr
374.5798740116257 yr
436.78734548626034 yr
414.5809738299951 yr
380.24350245443225 yr
355.78172017799204 yr
424.68551500149505 yr
422.16946643772656 yr
384.09169948200656 yrAnd we take a horizontal slice at about 2000m.
plothorizontalslice(age_in_yrs, grd, depth=2000u"m", color=:magma)
Or look at the horiontal mean
plothorizontalmean(age_in_yrs, grd)
That's it for this tutorial... Good job!
This page was generated using Literate.jl.
- 1DeVries, T., & Holzer, M. (2019). Radiocarbon and helium isotope constraints on deep ocean ventilation and mantle‐³He sources. Journal of Geophysical Research: Oceans, 124, 3036–3057. doi:10.1029/2018JC014716