Public Documentation

Documentation for FourierFlows.jl's public interface.

See the Internals section of the manual for internal package docs covering all submodules.

FourierFlows

FourierFlows.FourierFlowsModule
source

Problem

Domain

FourierFlows.OneDGridType
OneDGrid(dev::Device = CPU();
         nx, Lx,
         x0 = -Lx/2, nthreads = Sys.CPU_THREADS, effort = FFTW.MEASURE, 
         T = Float64, aliased_fraction = 1/3)

Construct a OneDGrid on device; by default on CPU(). Grid size is Lx, resolution is nx, and leftmost position is x0. FFT plans are generated for nthreads CPUs using FFTW flag effort. The float type is T. The aliased_fraction keyword determines the highest wavenubers that are being zero-ed out by dealias! function; 1/3 is the nominal value for quadratic nonlinearities.

source
FourierFlows.OneDGridType
struct OneDGrid{T<:AbstractFloat, A, Tx, Tfft, Trfft, Talias} <: AbstractGrid{T, A, Talias, D}

A one-dimensional grid.

  • device::Any: device which the grid lives on

  • nx::Int64: number of points in $x$

  • nk::Int64: number of wavenumbers in $x$

  • nkr::Int64: number of positive wavenumbers in $x$ (real Fourier transforms)

  • dx::AbstractFloat: grid spacing in $x$

  • Lx::AbstractFloat: domain extent in $x$

  • x::Any: range with $x$-grid-points

  • k::Any: array with $x$-wavenumbers

  • kr::Any: array with positive $x$-wavenumbers (real Fourier transforms)

  • invksq::Any: array with inverse squared $k$-wavenumbers, $1 / k²$

  • invkrsq::Any: array with inverse squared $kᵣ$-wavenumbers, $1 / kᵣ²$

  • fftplan::Any: the FFT plan for complex-valued fields

  • rfftplan::Any: the FFT plan for real-valued fields

  • aliased_fraction::AbstractFloat: the fraction of wavenumbers that are aliased (e.g., 1/3 for quadradic nonlinearities)

  • kalias::Any: range of the indices of aliased $x$-wavenumbers

  • kralias::Any: range of the indices of aliased positive $x$-wavenumbers (real Fourier transforms)

source
FourierFlows.ThreeDGridType
ThreeDGrid(dev::Device=CPU(); nx, Lx, ny=nx, Ly=Lx, nz=nx, Lz=Lx,
           x0=-Lx/2, y0=-Ly/2, z0=-Lz/2,
           nthreads=Sys.CPU_THREADS, effort=FFTW.MEASURE, T=Float64,
           aliased_fraction=1/3)

Construct a ThreeDGrid on device; by default on CPU(). Grid size is Lx, Ly, Lz, resolution is nx, ny, nz and leftmost positions are x0, y0, z0. FFT plans are generated for nthreads CPUs using FFTW flag effort. The float type is T. The aliased_fraction keyword determines the highest wavenubers that are being zero-ed out by dealias! function; 1/3 is the nominal value for quadratic nonlinearities.

source
FourierFlows.ThreeDGridType
struct ThreeDGrid{T<:AbstractFloat, Tk, Tx, Tfft, Trfft, Talias} <: AbstractGrid{T, Tk, Talias}

A three-dimensional grid.

  • device::Any: device which the grid lives on

  • nx::Int64: number of points in $x$

  • ny::Int64: number of points in $y$

  • nz::Int64: number of points in $z$

  • nk::Int64: number of wavenumbers in $x$

  • nl::Int64: number of wavenumbers in $y$

  • nm::Int64: number of wavenumbers in $z$

  • nkr::Int64: number of positive wavenumers in $x$ (real Fourier transforms)

  • dx::AbstractFloat: grid spacing in $x$

  • dy::AbstractFloat: grid spacing in $y$

  • dz::AbstractFloat: grid spacing in $z$

  • Lx::AbstractFloat: domain extent in $x$

  • Ly::AbstractFloat: domain extent in $y$

  • Lz::AbstractFloat: domain extent in $z$

  • x::Any: range with $x$-grid-points

  • y::Any: range with $y$-grid-points

  • z::Any: range with $z$-grid-points

  • k::Any: array with $x$-wavenumbers

  • l::Any: array with $y$-wavenumbers

  • m::Any: array with $z$-wavenumbers

  • kr::Any: array with positive $x$-wavenumbers (real Fourier transforms)

  • Ksq::Any: array with squared total wavenumbers, $k² + l² + m²$

  • invKsq::Any: array with inverse squared total wavenumbers, $1 / (k² + l² + m²)$

  • Krsq::Any: array with squared total wavenumbers for real Fourier transforms, $kᵣ² + l² + m²$

  • invKrsq::Any: array with inverse squared total wavenumbers for real Fourier transforms, $1 / (kᵣ² + l² + m²)$

  • fftplan::Any: the FFT plan for complex-valued fields

  • rfftplan::Any: the FFT plan for real-valued fields

  • aliased_fraction::AbstractFloat: the fraction of wavenumbers that are aliased (e.g., 1/3 for quadradic nonlinearities)

  • kalias::Any: range of the indices of aliased $x$-wavenumbers

  • kralias::Any: range of the indices of aliased positive $x$-wavenumbers (real Fourier transforms)

  • lalias::Any: range of the indices of aliased $y$-wavenumbers

  • malias::Any: range of the indices of aliased $m$-wavenumbers

source
FourierFlows.TwoDGridType
TwoDGrid(dev::Device=CPU(); nx, Lx, ny=nx, Ly=Lx,
         x0=-Lx/2, y0=-Ly/2, nthreads=Sys.CPU_THREADS, effort=FFTW.MEASURE,
         T=Float64, aliased_fraction=1/3)

Construct a TwoDGrid on device; by default on CPU(). Grid size is Lx, Ly, resolution is nx, ny, and leftmost positions are x0, y0. FFT plans are generated for nthreads CPUs using FFTW flag effort. The float type is T. The aliased_fraction keyword determines the highest wavenubers that are being zero-ed out by dealias! function; 1/3 is the nominal value for quadratic nonlinearities.

source
FourierFlows.TwoDGridType
struct TwoDGrid{T<:AbstractFloat, A, Tx, Tfft, Trfft, Talias, D} <: AbstractGrid{T, Tk, Talias, D}

A two-dimensional grid.

  • device::Any: device which the grid lives on

  • nx::Int64: number of points in $x$

  • ny::Int64: number of points in $y$

  • nk::Int64: number of wavenumbers in $x$

  • nl::Int64: number of wavenumbers in $y$

  • nkr::Int64: number of positive wavenumers in $x$ (real Fourier transforms)

  • dx::AbstractFloat: grid spacing in $x$

  • dy::AbstractFloat: grid spacing in $y$

  • Lx::AbstractFloat: domain extent in $x$

  • Ly::AbstractFloat: domain extent in $y$

  • x::Any: range with $x$-grid-points

  • y::Any: range with $y$-grid-points

  • k::Any: array with $x$-wavenumbers

  • l::Any: array with $y$-wavenumbers

  • kr::Any: array with positive $x$-wavenumbers (real Fourier transforms)

  • Ksq::Any: array with squared total wavenumbers, $k² + l²$

  • invKsq::Any: array with inverse squared total wavenumbers, $1 / (k² + l²)$

  • Krsq::Any: array with squared total wavenumbers for real Fourier transforms, $kᵣ² + l²$

  • invKrsq::Any: array with inverse squared total wavenumbers for real Fourier transforms, $1 / (kᵣ² + l²)$

  • fftplan::Any: the FFT plan for complex-valued fields

  • rfftplan::Any: the FFT plan for real-valued fields

  • aliased_fraction::AbstractFloat: the fraction of wavenumbers that are aliased (e.g., 1/3 for quadradic nonlinearities)

  • kalias::Any: range of the indices of aliased $x$-wavenumbers

  • kralias::Any: range of the indices of aliased positive $x$-wavenumbers (real Fourier transforms)

  • lalias::Any: range of the indices of aliased $y$-wavenumbers

source
FourierFlows.gridpointsMethod
gridpoints(grid::OneDDGrid)
gridpoints(grid::TwoDGrid)
gridpoints(grid::ThreeDGrid)

Return the collocation points of the grid in 1D (X), 2D (X, Y) or 3D arrays (X, Y, Z).

source

Utilities

FourierFlows.device_arrayMethod
device_array(grid::AbstractGrid)

Return the proper array type according to the grid's device, i.e., Array for CPU and CuArray for GPU.

source
FourierFlows.device_arrayMethod
device_array(device::Device)
device_array(device::Device, T, dim)

Return the proper array type according to the device, i.e., Array for CPU and CuArray for GPU.

source
FourierFlows.innereltypeMethod
innereltype(x)

Recursively determine the 'innermost' type in by the collection x (which may be, for example, a collection of a collection).

source
FourierFlows.superzerosMethod
superzeros(T, A)

Returns an array like A, but full of zeros. If innereltype(A) can be promoted to T, then the innermost elements of the array will have type T.

source
FourierFlows.@devzerosMacro
@devzeros dev T dims a b c...

Create arrays of all zeros with element type T, size dims, and global names a, b, c (for example) on device dev.

source
FourierFlows.@superzerosMacro
@superzeros T a b c d...
@superzeros T dims b c d...

Generate arrays b, c, d... with the super-dimensions of a and innereltype T.

source
FourierFlows.@zerosMacro
@zeros T dims a b c...

Create arrays of all zeros with element type T, size dims, and global names a, b, c (for example). An arbitrary number of arrays may be created.

source

Diagnostics

FourierFlows.DiagnosticType
mutable struct Diagnostic{T, N} <: AbstractDiagnostic

A diagnostic that includes N elements of type T.

  • calc::Function: function that returns the diagnostic via calc(prob)

  • prob::FourierFlows.Problem: the relevant problem for this diagnostic

  • data::Vector: vector where the diagnostic time-series is saved

  • t::Vector{Float64}: vector with the times for which the diagnostic was saved

  • steps::Vector{Int64}: vector with the problem's step for which the diagnostic was saved

  • freq::Int64: integer denoting how often (every how many problem.steps) to save the diagnostic

  • i::Int64: integer denoting how many times the diagnostic.data was updated

source
FourierFlows.DiagnosticMethod
Diagnostic(calc, prob; freq=1, nsteps=100, ndata=ceil(Int, (nsteps+1)/freq))

Construct a diagnostic that stores the result of calc(prob) with frequency freq.

Keywords

  • freq: Diagnostic is saved every freq steps.
  • nsteps: The total number of steps in problem.
  • ndata: The number of diagnostics to be saved.
source

Output

FourierFlows.OutputType
struct Output

The composite type for output.

  • prob::FourierFlows.Problem: the relevant problem for the output

  • path::String: the path for the output file

  • fields::Dict{Symbol, Function}: the fields to be saved; the relevant problem for this diagnostic

source
FourierFlows.OutputMethod
Output(prob, filename, fieldtuples...)

Define output for prob with fields and functions that calculate the output in the list of tuples fieldtuples = (fldname, func)....

source
FourierFlows.savediagnosticMethod
savediagnostic(diagnostic, diagname, filename)

Save diagnostic to filename under name diagname. Only the computed diagnostic is saved, that is, everything up to diagnostic's iteration diagnostic.i.

source

Timesteppers

FourierFlows.AB3TimeStepperType
AB3TimeStepper(equation::Equation, dev::Device=CPU())

Construct a 3rd order Adams-Bashforth timestepper for equation on device dev.

source
FourierFlows.AB3TimeStepperType
struct AB3TimeStepper{T} <: AbstractTimeStepper{T}

A 3rd-order Adams-Bashforth timestepper for time-stepping ∂u/∂t = RHS(u, t) via:

uⁿ⁺¹ = uⁿ + dt/12 * (23 * RHS(uⁿ, tⁿ) - 16 * RHS(uⁿ⁻¹, tⁿ⁻¹) + 5 * RHS(uⁿ⁻², tⁿ⁻²))

Adams-Bashforth is a multistep method, i.e., it not only requires information from the n-th time-step (uⁿ) but also from the previous two timesteps (uⁿ⁻¹ and uⁿ⁻²). For the first two timesteps, it falls back to a forward Euler timestepping scheme:

uⁿ⁺¹ = uⁿ + dt * RHS(uⁿ, tⁿ)
source
FourierFlows.ETDRK4TimeStepperType
struct ETDRK4TimeStepper{T,TL} <: AbstractTimeStepper{T}

A 4th-order exponential-time-differencing Runge-Kutta timestepper for time-stepping ∂u/∂t = L * u + N(u). The scheme treats the linear term L exactly while for the nonlinear terms N(u) it uses a 4th-order Runge-Kutta scheme (RK4TimeStepper). That is,

uⁿ⁺¹ = exp(L * dt) * uⁿ + RK4(N(uⁿ))

For more info refer to [2].

source
FourierFlows.ETDRK4TimeStepperType
ETDRK4TimeStepper(equation::Equation, dt, dev::Device=CPU())

Construct a 4th-order exponential-time-differencing Runge-Kutta timestepper with timestep dt for equation on device dev.

source
FourierFlows.FilteredAB3TimeStepperType
FilteredAB3TimeStepper(equation::Equation, dev::Device=CPU(); filterkwargs...)

Construct a 3rd order Adams-Bashforth timestepper with spectral filtering for equation on device dev.

source
FourierFlows.FilteredETDRK4TimeStepperType
FilteredETDRK4TimeStepper(equation, dt; filterkwargs...)

Construct a 4th-order exponential-time-differencing Runge-Kutta timestepper with timestep dt and spectral filtering for equation on device dev.

source
FourierFlows.FilteredLSRK54TimeStepperType
FilteredRK4TimeStepper(equation::Equation, dev::Device=CPU(); filterkwargs...)

Construct a 4th-order 5-stages 2-storage Runge-Kutta timestepper with spectral filtering for equation on device dev.

source
FourierFlows.FilteredRK4TimeStepperType
FilteredRK4TimeStepper(equation::Equation, dev::Device=CPU(); filterkwargs...)

Construct a 4th-order Runge-Kutta timestepper with spectral filtering for equation on device dev.

source
FourierFlows.ForwardEulerTimeStepperType
struct ForwardEulerTimeStepper{T} <: AbstractTimeStepper{T}

A forward Euler timestepper for time-stepping ∂u/∂t = RHS(u, t) via:

uⁿ⁺¹ = uⁿ + dt * RHS(uⁿ, tⁿ)
source
FourierFlows.LSRK54TimeStepperType
LSRK54TimeStepper(equation::Equation, dev::Device=CPU())

Construct a 4th-order 5-stages low-storage Runge-Kutta timestepper for equation on device dev.

source
FourierFlows.LSRK54TimeStepperType
struct LSRK54TimeStepper{T} <: AbstractTimeStepper{T}

A 4th-order 5-stages 2-storage Runge-Kutta timestepper for time-stepping ∂u/∂t = RHS(u, t) via:

S² = 0

for i = 1:5
  S² = Aᵢ * S² + dt * RHS(uⁿ, t₀ + Cᵢ * dt)
  uⁿ += Bᵢ * S²
end

uⁿ⁺¹ = uⁿ

where Aᵢ, Bᵢ, and Cᵢ are the $A$, $B$, and $C$ coefficients from the LSRK table at the $i$-th stage. For details, refer to [6].

Usage

The LSRK54TimeStepper is slower than the RK4TimeStepper but has less memory footprint; half compared to RK4TimeStepper.

If your simulation is bound by performance then use RK4TimeStepper; if your simulation is bound by memory then consider using LSRK54TimeStepper.

source
FourierFlows.RK4TimeStepperType
RK4TimeStepper(equation::Equation, dev::Device=CPU())

Construct a 4th-order Runge-Kutta timestepper for equation on device dev.

source
FourierFlows.RK4TimeStepperType
struct RK4TimeStepper{T} <: AbstractTimeStepper{T}

A 4th-order Runge-Kutta timestepper for time-stepping ∂u/∂t = RHS(u, t) via:

uⁿ⁺¹ = uⁿ + dt/6 * (k₁ + 2 * k₂ + 2 * k₃ + k₄)

where

k₁ = RHS(uⁿ, tⁿ)
k₂ = RHS(uⁿ + k₁ * dt/2, tⁿ + dt/2)
k₃ = RHS(uⁿ + k₂ * dt/2, tⁿ + dt/2)
k₄ = RHS(uⁿ + k₃ * dt, tⁿ + dt)
Usage

If your simulation is limited by memory then consider switching to LSRK54TimeStepper. The LSRK54TimeStepper timestepper has about half the memory footprint compared to the RK4TimeStepper with a about 25%-30% performance trade off.

source
FourierFlows.TimeStepperFunction
TimeStepper(stepper, equation, dt=nothing, dev=CPU(); kw...)

Instantiate the Timestepper for equation with timestep dt and on the device. The keyword arguments, kw, are passed to the timestepper constructor.

source
FourierFlows.stepforward!Method
stepforward!(prob::Problem, diags, nsteps::Int)

Step forward prob for nsteps, incrementing diags along the way. diags may be a single Diagnostic or a Vector of Diagnostics.

source

Diffusion Testbed Module

FourierFlows.Diffusion.ProblemFunction
Problem(dev::Device = CPU();
                 nx = 128,
                 Lx = 2π,
                  κ = 0,
                 dt = 0.01,
            stepper = "RK4",           
   aliased_fraction = 0,
                  T = Float64)

Construct a constant diffusivity problem.

Keyword arguments

  • dev: (required) CPU() or GPU(); computer architecture used to time-step problem.
  • nx: Number of grid points in $x$-domain.
  • Lx: Extent of the $x$-domain.
  • κ: Diffusivity coefficient.
  • dt: Time-step.
  • stepper: The extent of the $y$-domain.
  • aliased_fraction: the fraction of high-wavenubers that are zero-ed out by dealias!().
  • T: Float64 or Float32; floating point type used for problem data.
source