Diagnostics

VortexPasta.DiagnosticsModule
Diagnostics

Contains tools for computing different diagnostics (total energy, energy spectra, ...) from simulation data.

source

Kinetic energy

VortexPasta.Diagnostics.kinetic_energyFunction
kinetic_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.

source
VortexPasta.Diagnostics.kinetic_energy_from_streamfunctionFunction
kinetic_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

  • ψs: streamfunction values at filament nodes;

  • fs: vortex filament locations;

  • Γ::Real: quantum of circulation.

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. If nothing, 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.

source
VortexPasta.Diagnostics.kinetic_energy_nonperiodicFunction
kinetic_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 conservation

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.

source

Helicity

VortexPasta.Diagnostics.helicityFunction
helicity(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.

Arguments

  • fs: list of vortex filaments (or a single AbstractFilament);
  • vs: list of velocities on filament nodes (or a single vector of velocities on a filament);
  • Γ: vortex circulation.

See e.g. kinetic_energy_from_streamfunction for the meaning of the optional quad keyword argument.

source

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_rateFunction
stretching_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. If nothing, only values at nodes are used (cheaper).
source

Vortex impulse

VortexPasta.Diagnostics.vortex_impulseFunction
vortex_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).

source

Energy spectrum

VortexPasta.Diagnostics.energy_spectrum!Function
energy_spectrum!(Ek::AbstractVector, ks::AbstractVector, cache; unfilter = true)

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 a LongRangeCache);
  • a VortexFilamentSolver from the Timestepping module (which contains a BiotSavartCache).

The energy spectrum is computed from a recent Biot–Savart calculation using fast Ewald summation. More precisely, it is computed from the long-range velocity field in Fourier space. The LongRangeCache associated to the calculation is expected to currently contain this field.

In its most usual state, a LongRangeCache contains the long-range velocity field in the Ewald method, which is a Gaussian-filtered field (see e.g. BiotSavart.to_smoothed_velocity!). By default this function undoes the Gaussian filter, so that the returned kinetic energy spectrum is that of the unsmoothed velocity (which is singular at vortex positions, so it presents a slow decay in wavenumber space). One can pass unfilter = false to return the spectrum associated to the smoothed field.

The cache can also contain an unsmoothed 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 (unsmoothed) velocity field. Currently, the unfilter argument is ignored in this case.

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_energy_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.

source
VortexPasta.Diagnostics.energy_spectrumFunction
energy_spectrum(iter::VortexFilamentSolver; unfilter = true) -> (ks, Ek)
energy_spectrum(cache; unfilter = true) -> (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.

source
VortexPasta.Diagnostics.init_energy_spectrumFunction
Diagnostics.init_energy_spectrum(cache) -> (ks, Ek)

Initialise fields for storing an energy spectrum.

Returns a wavenumber vector ks and an uninitialised energy spectrum Ek with the right dimensions, which can be then passed to energy_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.

source