Velocity Interpolation

The dxdt! etc functions compute the tracked individual velocity.

IndividualDisplacements.dxdt!Function
dxdt!(du,u,p::𝐹_MeshArray3D,tim)

Interpolate velocity from gridded fields (3D; with halos) to position u (x,y,z,fIndex) to compute the derivative of position v time du_dt.

Extended help

using IndividualDisplacements
u,v,w,pos,func=vortex_flow_field(format=:MeshArray)
𝐹=FlowFields(u,u,v,v,0*w,1*w,[0,3*pi],func)
𝐼=Individuals(𝐹,pos...)
∫!(𝐼)

ref=[6.16, 7.23, 1.29, 1.0] # hide
prod(isapprox.(𝐼.📌,ref,atol=1.0)) # hide
source
dxdt!(du,u,p::𝐹_MeshArray2D,tim)

Interpolate velocity from gridded fields (2D; with halos) to position u (x,y,fIndex) to compute the derivative of position v time du_dt.

Extended help

using IndividualDisplacements
u,v,w,pos,func=random_flow_field(format=:MeshArray);
𝐹=FlowFields(u,u,v,v,[0,1.0],func);
𝐼=Individuals(𝐹,pos...);
∫!(𝐼);

isa(𝐼.📌,Vector)
source
dxdt!(du,u,𝑃::𝐹_Array3D,tim)

Interpolate velocity from gridded fields (3D; NO halos) to position u (x,y,z) to compute the derivative of position v time du_dt.

Extended help

using IndividualDisplacements
u,v,w,pos=vortex_flow_field(format=:Array)
𝐹=FlowFields(u,u,v,v,0*w,1*w,[0,3*pi])
𝐼=Individuals(𝐹,pos...)
∫!(𝐼)

ref=[9.35, 7.93, 1.28] # hide
prod(isapprox.(𝐼.📌,ref,atol=1.0)) # hide
source
dxdt!(du,u,𝑃::𝐹_Array2D,tim)

Interpolate velocity from gridded fields (2D; NO halos) to position u (x,y) to compute the derivative of position v time du_dt.

Extended help

using IndividualDisplacements
u,v,w,pos=random_flow_field(format=:Array);
𝐹=FlowFields(u,u,v,v,[0,1.0]);
𝐼=Individuals(𝐹,pos...);
∫!(𝐼);

isa(𝐼.📌,Vector)
source
IndividualDisplacements.dxy_dt_replayFunction
dxy_dt_replay(du,u,p::DataFrame,t)

Interpolate velocity from MITgcm float_trajectories output and return position increment du.

Extended help

using IndividualDisplacements
p=dirname(pathof(IndividualDisplacements))
include(joinpath(p,"../examples/more/detailed_look.jl"))
prod(isapprox.(sol[:,end],ref[:,end],atol=1.0))
source
IndividualDisplacements.dxy_dt_CyclicArrayFunction
dxy_dt_CyclicArray(du,u,𝑃::NamedTuple,tim)

Nearest neighbor (?) velocity from gridded fields (2D; NO halos but not needed when CyclicArrays is used to extend valid indice ranges).

Extended help

notes: spatial interpolation & temporal interpolation are lacking

using IndividualDisplacements
p=dirname(pathof(IndividualDisplacements))
include(joinpath(p,"../examples/more/example_CyclicArray.jl"))
(x,y)=cyclicarray_example()

ref=[330.5 290.5]
prod(isapprox.([x[end] y[end]],ref,atol=1.0))
source

Setup And Postprocessing

Convenience functions to initialize a simulation and post-process the output are provided.

IndividualDisplacements.postprocess_xyFunction
postprocess_xy(sol,𝑃::FlowFields,𝐷::NamedTuple; id=missing, 𝑇=missing)

Copy sol to a DataFrame & map position to x,y coordinates, and define time axis for a simple doubly periodic domain

source
IndividualDisplacements.postprocess_MeshArrayFunction
postprocess_MeshArray(sol,𝑃::FlowFields,𝐷::NamedTuple; id=missing, 𝑇=missing)

Copy sol to a DataFrame & map position to lon,lat coordinates using "exchanged" 𝐷.XC, 𝐷.YC via add_lonlat!

source
IndividualDisplacements.add_lonlat!Function
add_lonlat!(df::DataFrame,XC,YC)

Add lon & lat to dataframe using "exchanged" XC, YC

source
add_lonlat!(df::DataFrame,XC,YC,func::Function)

Add lon & lat to dataframe using "exchanged" XC, YC after updating subdomain indices (via func) if needed (MeshArrays.locationisout)

source

Basic geography support:

Base.diffFunction
Base.diff(𝐼::Individuals)

Difference in grid unit coordinates (dx,dy) between final and initial positions.

source
IndividualDisplacements.stprojFunction
stproj(XC,YC;XC0=0.0,YC0=90.0)

Apply to XC,YC (longitude, latitude) the stereographic projection which would put XC0,YC0 (longitude, latitude) at x,y=0,0

source
IndividualDisplacements.stproj_invFunction
stproj_inv(xx,yy;XC0=0.0,YC0=90.0)

Apply to xx,yy (stereographic projection coordinates) the reverse of the stereographic projection which puts XC0,YC0 (longitude, latitude) at x,y=0,0

source
IndividualDisplacements.interp_to_lonlatFunction
interp_to_lonlat(X::MeshArray,Γ::NamedTuple,lon,lat)

Use MeshArrays.Interpolate() to interpolate to arbitrary positions (e.g., a regular grid for plotting).

Extended help

using IndividualDisplacements
import IndividualDisplacements: MeshArrays
γ=MeshArrays.GridSpec("LatLonCap",MeshArrays.GRID_LLC90)
Γ=MeshArrays.GridLoad(γ,option="full")

lon=[i for i=20.:20.0:380., j=-70.:10.0:70.]
lat=[j for i=20.:20.0:380., j=-70.:10.0:70.]
tmp1=interp_to_lonlat(Γ.Depth,Γ,lon,lat)

(f,i,j,w,_,_,_)=MeshArrays.InterpolationFactors(Γ,vec(lon),vec(lat))
IntFac=(lon=lon,lat=lat,f=f,i=i,j=j,w=w)
tmp1=interp_to_lonlat(Γ.Depth,IntFac)
    
prod(isapprox(maximum(tmp1),5896.,atol=1.0))

# output

true
source
interp_to_lonlat(X::MeshArray,IntFac::NamedTuple)

Use MeshArrays.Interpolate() to interpolate to arbitrary positions (e.g., a regular grid for plotting).

source
IndividualDisplacements.nearest_to_xyFunction
nearest_to_xy(α::MeshArray,x,y,f)

Value of α at eachindex of the grid cell center nearest to x,y on subdomain array / facet f

source
nearest_to_xy(α::Array,x,y)

Value of α at eachindex of the grid cell center nearest to x,y

source

Toy Problems

These are used to demonstrate and test the package functionalities:

IndividualDisplacements.random_flow_fieldFunction
random_flow_field(;component=:Rotational,np=12,nq=18,format=:Array)

Generate random flow fields on a grid of np x nq points for use in simple examples.

The :Rotational component option is most similar to what is done in the standard example.

The :Divergent component option generates a purely divergent flow field instead.

(U,V,Φ)=random_flow_field(component=:Rotational)
F=convert_to_FlowFields(U,V,10.0)
I=Individuals(F,x,y,fill(1,length(x)))
source

Read External Files

Trajectory simulated by the MITgcm or observed by the global drifter program can be read from file using, respectively MITgcmTools.read_flt (from MITgcmTools.jl) or OceanRobots.drifters_hourly_read (from OceanRobots.jl).