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.FourierFlows
— ModuleMain module for FourierFlows.jl
– an ecosystem for solving partial differential equations on periodic domains using Fourier-based pseudospectral methods."
Exports
AB3TimeStepper
AbstractGrid
AbstractParams
AbstractTimeStepper
AbstractVars
CPU
Device
Diagnostic
ETDRK4TimeStepper
FilteredAB3TimeStepper
FilteredETDRK4TimeStepper
FilteredForwardEulerTimeStepper
FilteredLSRK54TimeStepper
FilteredRK4TimeStepper
ForwardEulerTimeStepper
GPU
LSRK54TimeStepper
OneDGrid
Output
RK4TimeStepper
ThreeDGrid
TimeStepper
TwoDGrid
@createarrays
cxtype
dealias!
device_array
@devzeros
fft
fltype
gridpoints
groupsize
ifft
increment!
innereltype
irfft
resize!
rfft
savediagnostic
saveoutput
saveproblem
step_until!
stepforward!
supersize
@superzeros
superzeros
update!
@zeros
FourierFlows.AbstractGrid
— TypeAbstract supertype for grids.
FourierFlows.AbstractParams
— TypeAbstract supertype for parameters.
FourierFlows.AbstractTimeStepper
— TypeAbstract supertype for timesteppers.
FourierFlows.AbstractVars
— TypeAbstract supertype for variables.
FourierFlows.CPU
— TypeCPU device.
FourierFlows.Device
— TypeAbstract supertype for supported devices.
FourierFlows.GPU
— TypeGPU device.
Problem
Domain
FourierFlows.OneDGrid
— TypeOneDGrid(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 dev
ice; 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.
FourierFlows.OneDGrid
— Typestruct OneDGrid{T<:AbstractFloat, A, Tx, Tfft, Trfft, Talias} <: AbstractGrid{T, A, Talias, D}
A one-dimensional grid
.
device::Any
: device which the grid lives onnx::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-pointsk::Any
: array with $x$-wavenumberskr::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 fieldsrfftplan::Any
: the FFT plan for real-valued fieldsaliased_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$-wavenumberskralias::Any
: range of the indices of aliased positive $x$-wavenumbers (real Fourier transforms)
FourierFlows.ThreeDGrid
— TypeThreeDGrid(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 dev
ice; 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.
FourierFlows.ThreeDGrid
— Typestruct ThreeDGrid{T<:AbstractFloat, Tk, Tx, Tfft, Trfft, Talias} <: AbstractGrid{T, Tk, Talias}
A three-dimensional grid
.
device::Any
: device which the grid lives onnx::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-pointsy::Any
: range with $y$-grid-pointsz::Any
: range with $z$-grid-pointsk::Any
: array with $x$-wavenumbersl::Any
: array with $y$-wavenumbersm::Any
: array with $z$-wavenumberskr::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 fieldsrfftplan::Any
: the FFT plan for real-valued fieldsaliased_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$-wavenumberskralias::Any
: range of the indices of aliased positive $x$-wavenumbers (real Fourier transforms)lalias::Any
: range of the indices of aliased $y$-wavenumbersmalias::Any
: range of the indices of aliased $m$-wavenumbers
FourierFlows.TwoDGrid
— TypeTwoDGrid(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 dev
ice; 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.
FourierFlows.TwoDGrid
— Typestruct 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 onnx::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-pointsy::Any
: range with $y$-grid-pointsk::Any
: array with $x$-wavenumbersl::Any
: array with $y$-wavenumberskr::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 fieldsrfftplan::Any
: the FFT plan for real-valued fieldsaliased_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$-wavenumberskralias::Any
: range of the indices of aliased positive $x$-wavenumbers (real Fourier transforms)lalias::Any
: range of the indices of aliased $y$-wavenumbers
FourierFlows.dealias!
— Methoddealias!(fh, grid)
Dealias array fh
on the grid
based on the grids
's aliased_fraction
.
FourierFlows.gridpoints
— Methodgridpoints(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
).
Utilities
FourierFlows.cxtype
— Methodcxtype(T)
Returns T
when T
is Complex
, or Complex{T}
when T
is Real
.
FourierFlows.device_array
— Methoddevice_array(grid::AbstractGrid)
Return the proper array type according to the grid
's device
, i.e., Array
for CPU and CuArray
for GPU.
FourierFlows.device_array
— Methoddevice_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.
FourierFlows.fltype
— Methodfltype(T)
Returns T
when T<:AbstractFloat
or Tf
when T<:Complex{Tf}
.
FourierFlows.innereltype
— Methodinnereltype(x)
Recursively determine the 'innermost' type in by the collection x
(which may be, for example, a collection of a collection).
FourierFlows.superzeros
— Methodsuperzeros(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
.
FourierFlows.@devzeros
— Macro@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
.
FourierFlows.@superzeros
— Macro@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
.
FourierFlows.@zeros
— Macro@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.
Diagnostics
FourierFlows.Diagnostic
— Typemutable struct Diagnostic{T, N} <: AbstractDiagnostic
A diagnostic that includes N
elements of type T
.
calc::Function
: function that returns the diagnostic viacalc(prob)
prob::FourierFlows.Problem
: the relevant problem for this diagnosticdata::Vector
: vector where the diagnostic time-series is savedt::Vector{Float64}
: vector with the times for which the diagnostic was savedsteps::Vector{Int64}
: vector with the problem'sstep
for which the diagnostic was savedfreq::Int64
: integer denoting how often (every how manyproblem.step
s) to save the diagnostici::Int64
: integer denoting how many times thediagnostic.data
was updated
FourierFlows.Diagnostic
— MethodDiagnostic(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 everyfreq
steps.nsteps
: The total number of steps in problem.ndata
: The number of diagnostics to be saved.
FourierFlows.increment!
— Methodincrement!(diag)
Increment the Diagnostic diag
, or an array of Diagnostics diags
.
FourierFlows.update!
— Methodupdate!(diag)
Update diag
.
Output
FourierFlows.Output
— Typestruct Output
The composite type for output.
prob::FourierFlows.Problem
: the relevant problem for the outputpath::String
: the path for the output filefields::Dict{Symbol, Function}
: the fields to be saved; the relevant problem for this diagnostic
FourierFlows.Output
— MethodOutput(prob, filename, fieldtuples...)
Define output for prob
with fields and functions that calculate the output in the list of tuples fieldtuples = (fldname, func)...
.
FourierFlows.savediagnostic
— Methodsavediagnostic(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
.
FourierFlows.saveoutput
— Methodsaveoutput(out)
Save the fields in out.fields
to out.path
.
FourierFlows.saveproblem
— Methodsaveproblem(prob, filename)
Save certain aspects of a problem prob
to filename
.
Timesteppers
FourierFlows.AB3TimeStepper
— TypeAB3TimeStepper(equation::Equation, dev::Device=CPU())
Construct a 3rd order Adams-Bashforth timestepper for equation
on device dev
.
FourierFlows.AB3TimeStepper
— Typestruct 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ⁿ)
FourierFlows.ETDRK4TimeStepper
— TypeETDRK4TimeStepper(equation::Equation, dt, dev::Device=CPU())
Construct a 4th-order exponential-time-differencing Runge-Kutta timestepper with timestep dt
for equation
on device dev
.
FourierFlows.ETDRK4TimeStepper
— Typestruct 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 Kassam and Trefethen (2005).
FourierFlows.FilteredAB3TimeStepper
— TypeFilteredAB3TimeStepper(equation::Equation, dev::Device=CPU(); filterkwargs...)
Construct a 3rd order Adams-Bashforth timestepper with spectral filtering for equation
on device dev
.
FourierFlows.FilteredAB3TimeStepper
— Typestruct FilteredAB3TimeStepper{T} <: AbstractTimeStepper{T}
A 3rd order Adams-Bashforth timestepper with spectral filtering. See AB3TimeStepper
.
FourierFlows.FilteredETDRK4TimeStepper
— TypeFilteredETDRK4TimeStepper(equation, dt; filterkwargs...)
Construct a 4th-order exponential-time-differencing Runge-Kutta timestepper with timestep dt
and spectral filtering for equation
on device dev
.
FourierFlows.FilteredETDRK4TimeStepper
— TypeFilteredETDRK4TimeStepper{T,TL,Tf} <: AbstractTimeStepper{T}
A 4th-order exponential-time-differencing Runge-Kutta timestepper with spectral filtering. See ETDRK4TimeStepper
.
FourierFlows.FilteredForwardEulerTimeStepper
— TypeFilteredForwardEulerTimeStepper(equation, dev; filterkwargs...)
Construct a Forward Euler timestepper with spectral filtering for equation
on device dev
.
FourierFlows.FilteredForwardEulerTimeStepper
— Typestruct FilteredForwardEulerTimeStepper{T,Tf} <: AbstractTimeStepper{T}
A forward Euler timestepper with spectral filtering. See ForwardEulerTimeStepper
.
FourierFlows.FilteredLSRK54TimeStepper
— TypeFilteredRK4TimeStepper(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
.
FourierFlows.FilteredLSRK54TimeStepper
— Typestruct FilteredLSRK54TimeStepper{T,V,Tf} <: AbstractTimeStepper{T}
A 4th-order 5-stages low-storage Runge-Kutta timestepper with spectral filtering. See LSRK54TimeStepper
.
FourierFlows.FilteredRK4TimeStepper
— TypeFilteredRK4TimeStepper(equation::Equation, dev::Device=CPU(); filterkwargs...)
Construct a 4th-order Runge-Kutta timestepper with spectral filtering for equation
on device dev
.
FourierFlows.FilteredRK4TimeStepper
— Typestruct FilteredRK4TimeStepper{T,Tf} <: AbstractTimeStepper{T}
A 4th-order Runge-Kutta timestepper with spectral filtering. See RK4TimeStepper
.
FourierFlows.ForwardEulerTimeStepper
— TypeForwardEulerTimeStepper(equation::Equation, dev::Device=CPU())
Construct a forward Euler timestepper for equation
on device dev
.
FourierFlows.ForwardEulerTimeStepper
— Typestruct ForwardEulerTimeStepper{T} <: AbstractTimeStepper{T}
A forward Euler timestepper for time-stepping ∂u/∂t = RHS(u, t)
via:
uⁿ⁺¹ = uⁿ + dt * RHS(uⁿ, tⁿ)
FourierFlows.LSRK54TimeStepper
— TypeLSRK54TimeStepper(equation::Equation, dev::Device=CPU())
Construct a 4th-order 5-stages low-storage Runge-Kutta timestepper for equation
on device dev
.
FourierFlows.LSRK54TimeStepper
— Typestruct 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 Carpenter and Kennedy (1994).
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
.
FourierFlows.RK4TimeStepper
— TypeRK4TimeStepper(equation::Equation, dev::Device=CPU())
Construct a 4th-order Runge-Kutta timestepper for equation
on device dev
.
FourierFlows.RK4TimeStepper
— Typestruct 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)
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.
FourierFlows.TimeStepper
— FunctionTimeStepper(stepper, equation, dt=nothing, dev=CPU(); kw...)
Instantiate the Timestepper
for equation
with timestep dt
and on the dev
ice. The keyword arguments, kw
, are passed to the timestepper constructor.
FourierFlows.step_until!
— Methodstep_until!(prob, stop_time)
Step forward prob
until stop_time
.
We cannot use step_until!
with ETDRK4TimeStepper
nor FilteredETDRK4TimeStepper
.
See also: stepforward!
FourierFlows.stepforward!
— Methodstepforward!(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 Diagnostic
s.
FourierFlows.stepforward!
— Methodstepforward!(prob::Problem, nsteps::Int)
Step forward prob
for nsteps
.
FourierFlows.stepforward!
— Methodstepforward!(prob::Problem)
Step forward prob
one time step.
Diffusion Testbed Module
FourierFlows.Diffusion
— ModuleFourierFlows.Diffusion.Problem
— FunctionProblem(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()
orGPU()
; computer architecture used to time-stepproblem
.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 bydealias!()
.T
:Float64
orFloat32
; floating point type used forproblem
data.
FourierFlows.Diffusion.set_c!
— Methodset_c!(prob, c)
Set the solution sol
as the transform of c
and update vars
.
FourierFlows.Diffusion.updatevars!
— Methodupdatevars!(vars, grid, sol)
Update the variables in vars
on the grid
with the solution in sol
.