Diagnostics
VortexPasta.Diagnostics
— ModuleDiagnostics
Contains tools for computing different diagnostics (total energy, energy spectra, ...) from simulation data.
Kinetic energy
VortexPasta.Diagnostics.kinetic_energy
— Functionkinetic_energy(iter::VortexFilamentSolver; quad = nothing) -> Real
kinetic_energy(fs, ψs, Γ::Real, [Ls]; quad = nothing) -> Real
kinetic_energy(fs, ψs, p::ParamsBiotSavart; quad = nothing) -> Real
Compute kinetic energy of velocity field induced by a set of vortex filaments.
This function simply calls kinetic_energy_from_streamfunction
.
VortexPasta.Diagnostics.kinetic_energy_from_streamfunction
— Functionkinetic_energy_from_streamfunction(iter::VortexFilamentSolver; quad = nothing)
kinetic_energy_from_streamfunction(fs, ψs, Γ::Real, [Ls]; quad = nothing)
kinetic_energy_from_streamfunction(fs, ψs, p::ParamsBiotSavart; quad = nothing)
Compute kinetic energy per unit mass (units $L^2 T^{-2}$) from streamfunction values at filament nodes in a periodic domain.
The kinetic energy per unit mass of the velocity field induced by a set of vortex filaments can be obtained as:
\[E = \frac{Γ}{2V} ∮ \bm{ψ}(\bm{s}) ⋅ \mathrm{d}\bm{s}\]
where $Γ$ is the vortex circulation and $V$ is the volume of interest (e.g. the volume of a periodic cell in periodic domains).
Arguments
Mandatory arguments
fs
: vortex filament locations;ψs
: streamfunction values at filament nodes;Γ::Real
: quantum of circulation; orp::ParamsBiotSavart
: Biot–Savart parameters (seeParamsBiotSavart
).
Optional arguments
Ls::Tuple = (Lx, Ly, Lz)
: domain size in each direction. If not given, the domain volume $V$ is taken to be 1 (see Non-periodic domains below).
Optional keyword arguments
quad = nothing
: optional quadrature rule (e.g.quad = GaussLegendre(4)
) used to evaluate line integrals. Ifnothing
, only values at nodes are used (cheaper). Otherwise, if a quadrature rule is passed, interpolations are performed and extra allocations are needed.
Non-periodic domains
If the domain is not periodic, that is, if one or more values in Ls
is Infinity
, then the domain volume $V$ in the expression above will be set to 1. This corresponds to computing an energy per unit density (units $L^5 T^{-2}$) instead of per unit mass, meaning that one should multiply by the fluid density $ρ$ ($M L^{-3}$) to get an actual kinetic energy ($M L^2 T^{-2}$).
Note that in non-periodic domains one may also use the kinetic_energy_nonperiodic
function, which uses a different definition commonly used for open domains and which does not require streamfunction values (but only velocity values on filament nodes). However, energy computed using that definition may not be properly conserved when it should.
Therefore, it is recommended to always use kinetic_energy_from_streamfunction
, even in non-periodic domains.
VortexPasta.Diagnostics.kinetic_energy_nonperiodic
— Functionkinetic_energy_nonperiodic(iter::VortexFilamentSolver; quad = nothing) -> Real
kinetic_energy_nonperiodic(fs, vs, Γ::Real; quad = nothing) -> Real
kinetic_energy_nonperiodic(fs, vs, p::ParamsBiotSavart; quad = nothing) -> Real
Compute kinetic energy per unit density (units $L^5 T^{-2}$) from velocity values at filament nodes in an open (non-periodic) domain.
This function returns the kinetic energy over the infinite fluid volume, which only makes sense in an open domain such that the velocity tends to zero far from the vortices.
In an open domain, the kinetic energy of the velocity field induced by a set of vortex filaments can be obtained as:
\[E = ρ Γ ∮ \bm{v} ⋅ (\bm{s} × \mathrm{d}\bm{s})\]
where $Γ$ is the vortex circulation and $ρ$ is the fluid density. This definition assumes that the velocity field tends to zero far from the vortices. This is true in open domains, but not in periodic ones.
Energy computed using this definition may present small temporal fluctuations in cases where energy should be conserved. For this reason, it is recommended to always use kinetic_energy_from_streamfunction
, which displays proper energy conservation properties. In general, this definition will slightly overestimate the energy obtained from the streamfunction.
This function returns the energy per unit density $E / ρ$.
Mandatory arguments
vs
: velocity values at filament nodes;fs
: vortex filament locations;Γ::Real
: quantum of circulation.
Optional keyword arguments
See kinetic_energy_from_streamfunction
for details.
Periodic domains
In periodic domains this function will give wrong results, since there are boundary terms coming from integration by parts which are neglected in the above definition (assuming the velocity goes to zero far from the vortices, which is not the case in a periodic domain).
In this case, the kinetic_energy_from_streamfunction
function should be used instead.
Energy injection rate
VortexPasta.Diagnostics.energy_injection_rate
— Functionenergy_injection_rate(iter::VortexFilamentSolver, [vL]; quad = nothing) -> Real
energy_injection_rate(fs, vL, vs, p::ParamsBiotSavart; quad = nothing) -> Real
Compute energy injection rate from current filament velocities.
This is typically the energy injected by an advecting vortex velocity vL
, given the self-induced (Biot–Savart) velocity vs
. It can be negative if vL
actually dissipates energy. It does not include the contributions to the dissipation by vortex reconnections or numerical resolution effects.
Arguments
fs
: vortex filament locations;vL
: velocity of filament nodes (including forcing or external velocities);vs
: self-induced filament velocity (from Biot–Savart law);p::ParamsBiotSavart
: Biot–Savart parameters (seeParamsBiotSavart
);quad
: an optional quadrature rule (such asGaussLegendre
) for accurate line integrals.
Extended help
Definition
The energy injection rate (in $L^2 T^{-3}$ units) can be expressed as
\[ε_{\text{inj}} = \frac{Γ}{V} ∮ \left[ \bm{s}' × \bm{v}_{\text{s}} \right] ⋅ \bm{v}_{\text{L}} \, \mathrm{d}ξ,\]
where $\bm{v}_{\text{s}}$ is the self-induced vortex velocity (according to the Biot–Savart law) and $\bm{v}_{\text{L}}$ is the actual velocity used to advect the vortices. Typically, the latter can be written as
\[\bm{v}_{\text{L}} = \bm{v}_{\text{s}} + \bm{v}_{\text{ext}}\]
where $\bm{v}_{\text{ext}}$ is an externally applied velocity (for example, representing the interaction with a normal fluid). Note that $\bm{v}_{\text{s}}$ doesn't contribute to energy injection, so that energy is conserved (excluding other dissipative effects) in the absence of an external velocity.
Helicity
VortexPasta.Diagnostics.helicity
— Functionhelicity(iter::VortexFilamentSolver; quad = nothing) -> Real
helicity(fs, vs, Γ::Real; quad = nothing) -> Real
helicity(fs, vs, p::ParamsBiotSavart; quad = nothing) -> Real
Compute helicity of a vortex configuration.
The helicity is defined as:
\[\mathcal{H} = ∫ \bm{u}(\bm{x}) ⋅ \bm{ω}(\bm{x}) \, \mathrm{d}^3\bm{x} = Γ ∮ \bm{v}(\bm{s}) ⋅ \mathrm{d}\bm{s}\]
where $\bm{u}(\bm{x})$ and $\bm{ω}(\bm{x}) = \bm{∇} × \bm{u}(\bm{x})$ are the velocity and vorticity fields, and $\bm{v}(\bm{s})$ is the velocity of the set of vortex filaments.
Note that the returned helicity has units of a squared circulation ($Γ^2 ∼ L^4 T^{-2}$).
Arguments
fs
: list of vortex filaments (or a singleAbstractFilament
);vs
: list of velocities on filament nodes (or a single vector of velocities on a filament);p
: Biot–Savart parameters; orΓ
: vortex circulation.
See e.g. kinetic_energy_from_streamfunction
for the meaning of the optional quad
keyword argument.
Filament length
See filament_length
in the Filaments
module. For convenience, filament_length
is re-exported by Diagnostics
, meaning that one can do:
using VortexPasta.Diagnostics
filament_length(...)
without needing to import VortexPasta.Filaments
.
Stretching rate
VortexPasta.Diagnostics.stretching_rate
— Functionstretching_rate(iter::VortexFilamentSolver; quad = nothing) -> Real
stretching_rate(fs, vs; quad = nothing) -> Real
Compute stretching rate of one or more vortices.
The stretching rate has units $L T^{-1}$ and is given by:
\[\frac{\mathrm{d} \mathcal{L}}{\mathrm{d} t} = ∮ \frac{∂\bm{v}}{∂ξ} ⋅ \mathrm{d}\bm{s} = - ∮ \bm{v} ⋅ \bm{s}'' \, \mathrm{d}ξ\]
where $ξ$ is the arc length, $\bm{v}(ξ)$ the local filament velocity, $\bm{s}''(ξ)$ the local curvature vector, and $\mathcal{L}$ the instantaneous vortex length. The last equality is obtained using integration by parts.
In the implementation, the last expression is the one used to compute the stretching rate.
Mandatory arguments
vs
: velocity values at filament nodes;fs
: vortex filament locations.
Optional keyword arguments
quad = nothing
: optional quadrature rule (e.g.quad = GaussLegendre(4)
) used to evaluate line integrals. Ifnothing
, only values at nodes are used (cheaper).
Vortex impulse
VortexPasta.Diagnostics.vortex_impulse
— Functionvortex_impulse(iter::VortexFilamentSolver; quad = nothing) -> Vec3
vortex_impulse(f; quad = nothing) -> Vec3
Estimate normalised impulse of one or more vortex filaments.
The vortex impulse is defined as
\[\bm{p} = \frac{1}{2} ρ Γ ∮ \bm{s} × \mathrm{d}\bm{s}\]
where $ρ$ is the fluid density and $Γ$ the circulation about the vortex. Note that this function returns the normalised impulse $\bm{p} / ρΓ$. The returned impulse has units of $L^2$ (an area).
Note that, for a circular vortex ring of radius $R$, its impulse is $\bm{p} = ρΓA$ where $A = π R^2$ is the area enclosed by the ring (and the orientation is equal to its direction of propagation, i.e. normal to the plane where the ring lives).
Energy and helicity spectrum
VortexPasta.Diagnostics.init_spectrum
— FunctionDiagnostics.init_spectrum(iter::VortexFilamentSolver) -> (ks, Ek)
Diagnostics.init_spectrum(cache) -> (ks, Ek)
Initialise fields for storing an energy or helicity spectrum.
Returns a wavenumber vector ks
and an uninitialised spectrum Ek
with the right dimensions, which can be then passed to energy_spectrum!
or `helicity_spectrum!.
The returned arrays are always on the CPU, even when the cache
contains GPU data.
See energy_spectrum!
for details on the cache
argument.
VortexPasta.Diagnostics.energy_spectrum
— Functionenergy_spectrum(iter::VortexFilamentSolver) -> (ks, Ek)
energy_spectrum(cache) -> (ks, Ek)
Compute kinetic energy spectrum associated to vortex filament state.
Returns a tuple of vectors (ks, Ek)
where ks
contains the probed wavenumbers and Ek
the energy associated to each wavenumber.
See also energy_spectrum!
for a non-allocating variant and for more details.
VortexPasta.Diagnostics.energy_spectrum!
— Functionenergy_spectrum!(Ek::AbstractVector, ks::AbstractVector, cache)
Compute kinetic energy spectrum associated to vortex filament state.
Here cache
contains the results of long-range Biot–Savart computations. It can be either:
- a
LongRangeCache
; - a
BiotSavartCache
(which contains aLongRangeCache
); - a
VortexFilamentSolver
from theTimestepping
module (which contains aBiotSavartCache
).
The energy spectrum is computed from a recent Biot–Savart calculation using fast Ewald summation. More precisely, it is computed from the (Fourier-truncated) velocity field in Fourier space. The LongRangeCache
associated to the calculation is expected to currently contain this field.
The cache can also contain the vorticity field in Fourier space (the result obtained right after performing a NUFFT from the filament locations, see BiotSavart.compute_vorticity_fourier!
. In this case this function does the right thing and also computes the spectrum of the associated velocity field.
The vectors Ek
and ks
are expected to have the same length. Moreover, the vector of wavenumbers ks
should satisfy ks[begin] == 0
and have a constant step Δk = ks[i + 1] - ks[i]
. For convenience, the init_spectrum
function can be used to create these vectors.
See also energy_spectrum
for an allocating variant which doesn't need predefined Ek
and ks
vectors.
VortexPasta.Diagnostics.helicity_spectrum
— Functionhelicity_spectrum(iter::VortexFilamentSolver) -> (ks, Hk)
helicity_spectrum(cache) -> (ks, Hk)
Compute helicity spectrum associated to vortex filament state.
Returns a tuple of vectors (ks, Hk)
where ks
contains the probed wavenumbers and Hk
the helicity associated to each wavenumber.
See also helicity_spectrum!
for a non-allocating variant and for more details.
VortexPasta.Diagnostics.helicity_spectrum!
— Functionhelicity_spectrum!(Hk::AbstractVector, ks::AbstractVector, cache)
Compute helicity spectrum associated to vortex filament state.
For a single wave vector $\bm{k}$, the associated helicity is defined as
\[H(\bm{k}) ≡ \frac{V}{Δk} \bm{v}^*_{\bm{k}} ⋅ \bm{ω}_{\bm{k}} = 2\frac{V}{Δk} \bm{k} ⋅ \left[ \bm{v}_{\text{r}} × \bm{v}_{\text{i}} \right].\]
In the last expression, the subscripts $\text{r}$ and $\text{i}$ respectively denote real and imaginary parts. That result can be easily found by doing some algebra, after replacing $\bm{ω}_{\bm{k}} = \mathrm{i} \bm{k} × \bm{u}_{\bm{k}}$. It explicitly shows that $H(\bm{k})$ is purely real.
As usual, this function computes the sum of $H(\bm{k})$ over spherical shells of radius $k_n ≤ |\bm{k}| < k_n + Δk = k_{n + 1}$.
Defined as above, the total helicity in a (cubic) periodic domain of volume $V = L^3$ is $H = ∑_{\bm{k}} H(\bm{k}) Δk$ where $Δk = 2π/L$ is the resolution in wavenumber space.
See energy_spectrum!
for details on accepted cache
arguments.