Functions and types
Functions
ILMPostProcessing.adams_bashforth_2_backward
— Methodadams_bashforth_2_backward(initial_conditions, u, v, t0, t_start, dt, T)
Solve the initial value problem (IVP) using 2-step Adams-Bashforth method.
Integrate backward in time to compute backward FTLE fields.
Arguments
initial_conditions
: generated with function geninitcondsu
: array of x-velocity fieldsv
: array of y-velocity fieldst0
: initial timet_start
: start time ofu
andv
dt
: step size between consecutive velocity fieldsT
: length of integration time
ILMPostProcessing.adams_bashforth_2_forward
— Methodadams_bashforth_2_forward(initial_conditions, u, v, t0, t_start, dt, T)
Solve the initial value problem (IVP) using 2-step Adams-Bashforth method.
Integrate forward in time to compute forward FTLE fields.
Arguments
initial_conditions
: generated with function geninitcondsu
: array of x-velocity fieldsv
: array of y-velocity fieldst0
: initial timet_start
: start time ofu
andv
dt
: step size between consecutive velocity fieldsT
: length of integration time
ILMPostProcessing.compute_FTLE!
— Methodcompute_FTLE!(FTLE, nx, ny, T, final_positions, dx, dy)
Compute the FTLE
field given the final positions of initial points on a collocated grid.
The underlying math is detailed in: https://shaddenlab.berkeley.edu/uploads/LCS-tutorial/computation.html. For each grid point, first compute the gradient of the flow map using two point central differencing. Then, calculate the maximum eigenvalue of the 2 x 2 gradient matrix. Finally, compute the FTLE value using the eigenvalue.
Arguments
FTLE
: an empty 2D array (i.e.,FTLE = zeros(Float64, ny - 2, nx - 2))
,nx - 2
andny - 2
accounts for the boundary points in the central difference formulanx
: number of grid points in xny
: number of grid points in yT
: length of integration timefinal_positions
: solutions of the IVPdx
: spacing of initial grids in xdy
: spacing of initial grids in y
ILMPostProcessing.compute_FTLE!
— Methodcompute_FTLE!(FTLE::ScalarGridData, final_x::ScalarGridData, final_y::ScalarGridData, dx::Real, dy::Real, T::Real)
Compute the FTLE
field given the final positions of initial points on a collocated grid.
The underlying math is detailed in: https://shaddenlab.berkeley.edu/uploads/LCS-tutorial/computation.html. For each grid point, first compute the gradient of the flow map using two point central differencing. Then, calculate the maximum eigenvalue of the 2 x 2 gradient matrix. Finally, compute the FTLE value using the eigenvalue.
Arguments
FTLE
: will hold the FTLE field, inScalarGridData
typefinal_x
: deformed x positions, inScalarGridData
typefinal_y
: deformed y positions, inScalarGridData
typedx
: spacing of initial grids in xdy
: spacing of initial grids in yT
: length of integration time
ILMPostProcessing.compute_streakline
— Methodcompute_streakline(u,v,X₀::Vector,t[;τmin = t-3.0, dtstreak=0.01,dttraj=0.001]) -> Vector, Vector
Calculate a streakline at time t
for a velocity field u
and v
, based on an injection point X₀
. The end of the streakline is set by τmin
, the earliest time at which a particle passed through the injection point. It defaults to 3 time units before the current instant t
. The time step size dtstreak
sets the resolution of the streakline (i.e., how often the particles are sampled along the streakline). It returns arrays of the x and y coordinates of the streakline.
ILMPostProcessing.compute_streamline
— Methodcompute_streamline(u,v,X₀,srange::Tuple,t::Real,[;dt=0.01])
Calculate the streamline(s) passing through location(s) X₀
, which can be specified as either a single vector [x0,y0]
, a vector of vectors specifying x, y pairs, or a tuple of vectors or arrays specifying x and y positions, respectively, for a rake of streamlines. The arguments u
and v
are either interpolated velocity field components from a computational solution or are functions. If they are functions, then each of them should be of the form u(x,y,t)
and v(x,y,t)
; srange
is a tuple of the initial and final time of integration; t
is the current time at which the streamline is depicted; and Δs
is the time-like step size, which defaults to 0.001. The output is the solution structure for the OrdinaryDiffEq
package (or, for multiple points, a vector of such solution structures).
ILMPostProcessing.compute_trajectory
— Methodcompute_trajectory(v::VectorFieldSequence,X₀,Trange::Tuple[;dt=step(tr),alg=Euler()])
Calculate the trajectories of particles with initial location(s) X₀
. The argument v
contains a sequence of spatially-interpolated velocity fields and an associated time array. (One can also provide the vector of velocity fields and the time array as separate arguments).
X₀
can be specified as either a single vector [x0,y0]
, a vector of vectors specifying x, y pairs, or a tuple of vectors or arrays specifying x and y positions, respectively, for multiple tracer particles. Trange
is a tuple of the starting and ending integration times. The optional keyword arguments are dt
, the time step size (which defaults to the step size in tr
, but could be an integer multiple larger than 1). The output is the solution structure for the OrdinaryDiffEq
package.
ILMPostProcessing.compute_trajectory
— Methodcompute_trajectory(u,v,X₀,Trange::Tuple[;dt=0.01])
Calculate the trajectory of a tracer particle with initial location(s) X₀
, which can be specified as either a single vector [x0,y0]
, a vector of vectors specifying x, y pairs, or a tuple of vectors or arrays specifying x and y positions, respectively, for multiple tracer particles. The arguments u
and v
are either interpolated velocity field components from a computational solution or are functions. If they are functions, then each of them should be of the form u(x,y,t)
and v(x,y,t)
; Trange
is a tuple of the initial and final time of integration (and the final time can be earlier than the initial time if backward trajectories are desired); and dt
is the time step size, which defaults to 0.001. The output is the solution structure for the OrdinaryDiffEq
package (or, for multiple particles, a vector of such solution structures).
ILMPostProcessing.ddt
— Methodddt(u::AbstractVector,Δt[,mydiff=:backward_diff])
Calculate the time derivative of vector data u
, with time step size Δt
. The default method is backward differencing, but this can be changed to :forward_diff
or :central_diff
.
ILMPostProcessing.displacement_field
— Methoddisplacement_field(u::Function,v::Function,x0,y0,Trange::Tuple[;dt=:auto,alg=Euler()])
Calculate the displacement of particles initally at coordinates x0
and y0
over the range of times Trange = (ti,tf)
, using the velocity functions u
and v
. These function can either be autonomous (taking only x and y arguments) or non-autonomous, taking an additional time argument. The final time in Trange
can be earlier than the initial time if backward trajectories are desired. The optional keyword arguments are dt
, the time step size (which defaults to :auto
, for adaptive time marching, but could be specified to override this). The default time marching algorithm is Tsit5()
.
ILMPostProcessing.displacement_field
— Methoddisplacement_field(v::VectorFieldSequence,x0,y0,Trange::Tuple[;dt=step(tr),alg=Euler()])
Calculate the displacement of particles initally at coordinates x0
and y0
over the range of times Trange = (ti,tf)
, using the sequence of spatially-interpolated velocity fields in v
. (One can also provide the vector of velocity fields and the time array as separate arguments). The final time in Trange
can be earlier than the initial time if backward trajectories are desired. The optional keyword arguments are dt
, the time step size (which defaults to the step size in tr
, but could be an integer multiple larger than 1).
ILMPostProcessing.dmd
— Methoddmd(Xfull,r)
Compute the first r
DMD modes from the extended snapshot data in Xfull
. Both the original and shifted data are drawn from Xfull
, that is: X = Xfull[1:m-1]
and Xp = Xfull[2:m]
This returns the DMD modes and DMD eigenvalues in a DMDModes
structure, with fields modes
and evals
.
ILMPostProcessing.euler_backward
— Methodeuler_backward(initial_conditions, u, v, t0, t_start, dt, T)
Solve the initial value problem (IVP) using forward Euler method.
Integrate backward in time to compute backward FTLE fields. Note: not backward Euler method.
Arguments
initial_conditions
: generated with function geninitcondsu
: array of x-velocity fieldsv
: array of y-velocity fieldst0
: initial timet_start
: start time ofu
andv
dt
: step size between consecutive velocity fieldsT
: length of integration time
ILMPostProcessing.euler_forward
— Methodeuler_forward(initial_conditions, u, v, t0, t_start, dt, T)
Solve the initial value problem (IVP) using forward Euler method.
Integrate forward in time to compute forward FTLE fields.
Arguments
initial_conditions
: generated with function geninitcondsu
: array of x-velocity fieldsv
: array of y-velocity fieldst0
: initial timet_start
: start time ofu
andv
dt
: step size between consecutive velocity fieldsT
: length of integration time
ILMPostProcessing.field_along_trajectory
— Methodfieldalongtrajectory(f,traj::Trajectories,p::Integer[,deriv=0])
Evaluate field f
(given as grid data) along the trajectory number p
in the trajectories specified by traj
. The output is the history of f
along this trajectory. If f
is a vector field, then the component histories are output as a tuple. If deriv=1
, then it computes the time derivative of the field along the trajectory. The default is deriv=0
(no derivative).
ILMPostProcessing.gen_init_conds
— Methodgen_init_conds(X_MIN, X_MAX, Y_MIN, Y_MAX, nx, ny)
Generate a list of initial points (x, y).
These initial conditions represent a collocated grid with nx
points in (X_MIN
, X_MAX
) and ny
points in (Y_MIN
, Y_MAX
). The points are then flattened to a 1D array. The initial conditions could be used to compute FTLE or to visaulize trajectories.
ILMPostProcessing.make_interp_fields!
— Methodmake_interp_fields!(u, v, t_start, t_end, dt, velocity_fn, sol, sys, grid)
Generate an array of interpolatable velocity fields u and v using the solution of a ViscousFlow problem.
Each element of u
and v
is essentially a "Vectorfield!" used when solving an ODEProblem. Note: This function could be included in ViscousFlow.jl.
Arguments
u
: x-velocity fieldsv
: y-velocity fieldst_start
: start time of solt_end
: end time of soldt
: step size between consecutive velocity fieldsvelocity_fn
: function to compute velocity from solution, should be ViscousFlow.velocitysys
: viscous flow systemgrid
: physical grid
ILMPostProcessing.pod
— Methodpod(X::Vector{T}[; tolerance=0.99])
Calculate the POD modes and associated time-varying coefficients from an array of snapshot data X
. This X
plays the role of a snapshot matrix, whose columns are snapshots of the data. However, it is actually to be stored as a type Vector{T}
where T<:GridData
. It can be generated with a function like velocity(sol,sys,trange)
, where sol
is a ODESolution
, sys
is an ILMSystem
, and trange
is an array of times, e.g., trange=range(1,10,length=100)
. The number of POD modes retained in the decomposition is set by tolerance
: this specifies the fraction of the total energy to keep, and defaults to 99 percent.
The output of PODModes
is a structure with the following fields
Xmean
: temporal mean of the data. typeT
Xnorm
: originalX
vector with mean removed. Each element is of typeT
phi
: vector of POD modes. Each element is of typeT
a
: matrix of POD coefficients. Number of columns is same as number of entries inphi
. Columnk
constitutes the time-varying coefficient for modek
inphi
.lambda
: vector of modal energies, arranged in decreasing order, corresponding to the modes inphi
psi
: matrix of
ILMPostProcessing.pod
— Methodpod(X::Vector{T},r::Int)
Perform POD on snapshot data X
and truncate to r
modes
Types
ILMPostProcessing.ScalarFieldSequence
— TypeScalarFieldSequence(t::AbstractVector,s)
This type bundles a time vector t
with a vector s
of interpolatable scalar fields (i.e., each element is of type AbstractInterpolation
with two spatial coordinate arguments). It is used for plotting fields along trajectories.
ILMPostProcessing.Trajectories
— TypeTrajectories
Type returned by trajectory calculations, containing a set of one or more trajectories. For an instance traj
of this type, the trajectory time array can be returned with traj.t
. The number of trajectories is returned by traj.np
. Any trajectory contained in the set can be obtained with traj[p]
, where p
must be 0 < p <= traj.np
.
ILMPostProcessing.VectorFieldSequence
— TypeVectorFieldSequence(t::AbstractVector,v)
This type bundles a time vector t
with a vector v
of tuples of interpolatable fields (i.e., each member of the tuple is of type AbstractInterpolation
with two spatial coordinate arguments). It is used in trajectory computations and for plotting fields along trajectories.
Index
ILMPostProcessing.ScalarFieldSequence
ILMPostProcessing.Trajectories
ILMPostProcessing.VectorFieldSequence
ILMPostProcessing.adams_bashforth_2_backward
ILMPostProcessing.adams_bashforth_2_forward
ILMPostProcessing.compute_FTLE!
ILMPostProcessing.compute_FTLE!
ILMPostProcessing.compute_streakline
ILMPostProcessing.compute_streamline
ILMPostProcessing.compute_trajectory
ILMPostProcessing.compute_trajectory
ILMPostProcessing.ddt
ILMPostProcessing.displacement_field
ILMPostProcessing.displacement_field
ILMPostProcessing.dmd
ILMPostProcessing.euler_backward
ILMPostProcessing.euler_forward
ILMPostProcessing.field_along_trajectory
ILMPostProcessing.gen_init_conds
ILMPostProcessing.make_interp_fields!
ILMPostProcessing.pod
ILMPostProcessing.pod