Private types and functions
Documentation for FourierFlows.jl's internal interface.
FourierFlows
FourierFlows.AbstractDiagnostic — TypeAbstract supertype for diagnostics.
Problem
FourierFlows.Clock — Typemutable struct Clock{T<:AbstractFloat}Represents the clock of a problem.
dt::AbstractFloat: the time-stept::AbstractFloat: the timestep::Int64: the step number
FourierFlows.EmptyParams — Typestruct EmptyParams <: AbstractParamsA placeholder struct for parameters.
FourierFlows.EmptyVars — Typestruct EmptyVars <: AbstractVarsA placeholder struct for variables.
FourierFlows.Equation — Typestruct Equation{T, TL, G<:AbstractFloat}The equation to be solved ∂u/∂t = L * u + N(u). Array L includes the coefficients of the linear term L * u and calcN! is a function which computes the nonlinear term N(u). The struct also includes the problem's grid and the float type of the state vector (and consequently of N(u)).
L::Any: array with coefficient for the linear part of the equationcalcN!::Function: function that computes the nonlinear part of the equationgrid::AbstractGrid{G} where G<:AbstractFloat: the griddims::Tuple: the dimensions ofLT::Any: the float type for the state vector
FourierFlows.Equation — MethodEquation(L, calcN!, grid; dims=supersize(L), T=nothing)The equation constructor from the array L of the coefficients of the linear term, the function calcN! that computes the nonlinear term and the grid for the problem.
FourierFlows.Problem — Typestruct Problem{T, A<:AbstractArray, Tg<:AbstractFloat, TL}A problem that represents a partial differential equation.
sol::AbstractArray: the state vectorclock::FourierFlows.Clock: the problem'seqn::FourierFlows.Equation{T, TL, Tg} where {T, Tg<:AbstractFloat, TL}: the equationgrid::AbstractGrid{Tg} where Tg<:AbstractFloat: the gridvars::AbstractVars: the variablesparams::AbstractParams: the parameterstimestepper::AbstractTimeStepper{A} where A<:AbstractArray: the timestepper
FourierFlows.Problem — MethodProblem(eqn::Equation, stepper, dt, grid::AbstractGrid{T},
vars=EmptyVars, params=EmptyParams; stepperkwargs...) where TConstruct a Problem for equation eqn using the timestepper with timestep dt, on grid. The device is inferred from the grid. Optionally, use the keyword arguments to provide variables with vars and parameters with params. The stepperkwargs are passed on to the time-stepper constructor.
Domain
FourierFlows.getaliasedwavenumbers — Methodgetaliasedwavenumbers(nk, nkr, aliased_fraction)Return the top aliased_fraction highest wavenumbers, both for and real FFTs, kalias and kralias respectively. For example, aliased_fraction = 1/3 should return the indices of the top-most 1/6-th (in absolute value) for both positive and negative wavenumbers (i.e., 1/3 total) that should be set to zero after performing an FFT.
FourierFlows.makefilter — Methodmakefilter(K; order=4, innerK=2/3, outerK=1, tol=1e-15)Return a filter acting on the non-dimensional wavenumber K. For a one-dimensional grid, the non-dimensional wavenumber K is
K = k * dx / πand thus take values in $[-1, 1]$.
For K ≤ innerK the filter is inactive, i.e., equal to 1. For K > innerK, the filter decays exponentially to remove high-wavenumber content from the solution, i.e.,
filter(K) = exp(- decay * (K - innerK)^order)For a given order, the decay rate is determined so that the filter value at the outer wavenumber outerK is tol, where tol is a small number, close to machine precision.
decay = - log(tol) / (outerK - innerK)^orderUtilities
FourierFlows.jacobian — Methodjacobian(a, b, grid)Return the Jacobian of a and b on grid.
FourierFlows.jacobianh — Methodjacobianh(a, b, grid)Return the Fourier transform of the Jacobian of a and b on grid.
FourierFlows.on_grid — Methodon_grid(func, grid)Return an array, of the type compatible with the device that the grid lives on, that contains the values of function func evaluated on the grid.
FourierFlows.parsevalsum — Methodparsevalsum(uh, grid)Return the real part of the sum of uh on the grid. For example on a 2D grid, parsevalsum returns
\[ℜ [ \sum_{𝐤} û_{𝐤} L_x L_y ] , \]
where $û_{𝐤} =$ uh $/ (n_x e^{i 𝐤 ⋅ 𝐱₀})$. The elements of the vector $𝐱₀$ are the left-most position in each direction, e.g., for a 2D grid (grid.x[1], grid.y[1]).
When the input uh comes from a real-FFT transform, parsevalsum takes care to count the contribution from certain $k$-modes twice.
FourierFlows.parsevalsum2 — Methodparsevalsum2(uh, grid)Return the sum of |uh|² on the grid, which is equal to the domain integral of u². For example on a 2D grid, parsevalsum2 returns
\[\sum_{𝐤} |û_{𝐤}|² L_x L_y = \iint u² 𝖽x 𝖽y ,\]
where $û_{𝐤} =$ uh $/ (n_x e^{i 𝐤 ⋅ 𝐱₀})$. The elements of the vector $𝐱₀$ are the left-most position in each direction, e.g., for a 2D grid (grid.x[1], grid.y[1]).
When the input uh comes from a real-FFT transform, parsevalsum2 takes care to count the contribution from certain $k$-modes twice.
FourierFlows.radialspectrum — Methodradialspectrum(fh, grid; n=nothing, m=nothing, refinement=2)Return the radial spectrum of fh. fh lives on Cartesian wavenumber grid $(k, l)$. To compute the radial spectrum, we first interpolate $f̂(k, l)$ onto a radial wavenumber grid $(ρ, θ)$, where $ρ² = k² + l²$ and $θ = \arctan(l / k)$. Note here that $f̂ =$ fh $/ (n_x e^{i 𝐤 ⋅ 𝐱₀})$. The elements of the vector $𝐱₀$ are the left-most position in each direction, e.g., for a 2D grid (grid.x[1], grid.y[1]).
After interpolation, we integrate $f̂$ over angles $θ$ to get fρ,
\[ f̂_ρ = \int f̂(ρ, θ) ρ 𝖽ρ 𝖽θ .\]
The resolution (n, m) for the polar wavenumber grid is n = refinement * maximum(grid.nk, grid.nl), m = refinement * maximum(grid.nk, grid.nl), where refinement = 2 by default. If fh is in conjugate symmetric form then only the upper-half plane in $θ$ is represented on the polar grid.
Diagnostics
Base.getindex — Methode.g. plot(energydiag[:t], energydiag[:data]).
FourierFlows.extend! — Methodextend!(diag, n)Extend the data, time, and steps vectors of the diagnostic diag by n.
FourierFlows.extend! — Methodextend!(diag::Diagnostic{T,N})Double the extend of the data, time, and steps vectors of the diagnostic diag.
Output
FourierFlows.savefield — Methodsavefield(file, location, data)Saves a particular field's data to file.
FourierFlows.savefields — Methodsavefields(file, grid)Saves some parameters of problem's grid to file.
FourierFlows.uniquepath — Methoduniquepath(path)Return path with a number appended if isfile(path) == true. The number is incremented until path does not exist.
Timesteppers
FourierFlows.getetdcoeffs — Methodgetetdcoeffs(dt, L; ncirc=32, rcirc=1)Calculate the coefficients associated with the (diagonal) linear coefficient L for an ETDRK4 timestepper with timestep dt.
The calculation is done by integrating over a unit circle in the complex space. For more info refer to Kassam and Trefethen (2005).
Diffusion Testbed Module
FourierFlows.Diffusion.Params — Typestruct Params{T} <: AbstractParamsThe parameters for diffusion problem:
κ::Any: diffusivity coefficient
FourierFlows.Diffusion.Vars — Typestruct Vars{Aphys, Atrans} <: AbstractVarsThe variables for diffusion problem:
c: tracer concentrationcx: tracer concentration derivativech: Fourier transform of tracer concentrationcxh: Fourier transform of tracer concentration derivative
FourierFlows.Diffusion.Vars — MethodVars(grid)Return the variables vars for a constant diffusivity problem on grid.
FourierFlows.Diffusion.Equation — MethodEquation(grid, params)Return the equation for a constant diffusivity problem on grid with diffusivity found inside params.
FourierFlows.Diffusion.calcN! — MethodcalcN!(N, sol, t, clock, vars, params, grid)Calculate the nonlinear term for the 1D diffusion equation.