How to...
WaterWaves1D.jl
is meant to be versatile, and integrating new blocks to the package is easy. If you ever do so, please do not hesitate to contact the developers to either get help, or report on your advances.
build your model
Let us add the linear (Airy) water waves model whose equations are (using the notations introduced here)
\[ \left\{\begin{array}{l} ∂_tη+\frac{\tanh(\sqrt{μ} D)}{ν\sqrt{μ} D}∂_xv=0,\\[1ex] ∂_tv+∂_xη=0, \end{array}\right.\]
where we use the notation $F(D)$ for the action of pointwise multiplying by the function $F$ in the Fourier space.
In a dedicated file we write
export Airy
mutable struct Airy <: AbstractModel
label :: String
f! :: Function
mapto :: Function
mapfro :: Function
# We build our model here.
end
Our purpose is to provide
label
, a string used in subsequent informational messages, plots, etc.f!
, a function to be called in explicit time-integration solvers such asEuler
orRK4
(one may provide other functions to be used with other solvers such asEulerSymp
)mapto
, a function which from a couple(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executedmapfro
, a function which from raw data returns(η,v,x)
whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
.
To this aim, in place of the commented line, we write
function Airy(param::NamedTuple; # param is a NamedTuple containing all necessary parameters
label = "linear (Airy)" # using a keyword argument allows the user to supersede the default label.
)
# Set up
μ = param.μ
if !in(:ν,keys(param)) # set default ν if it is not provided
ν = min(1,1/√μ)
else
ν = param.ν
end
# Collocation points and Fourier modes
m = Mesh(param)
x, k = m.x, m.k
# Fourier multipliers
∂ₓ = 1im * k # Differentiation
∂ₓF₁ = 1im * tanh.(√μ*k)/(√μ*ν)
# Pre-allocation
fftη = zeros(Complex{Float64}, m.N)
fftv = zeros(Complex{Float64}, m.N)
# Evolution equations are ∂t U = f(U)
function f!(U)
fftη .= U[:,1]
fftv .= U[:,2]
U[:,1] .= -∂ₓF₁.*fftv
U[:,2] .= -∂ₓ.*fftη
end
# Build raw data from physical data (discrete Fourier transform)
function mapto(data::InitialData)
U = [fft(data.η(x)) fft(data.v(x))]
end
# Return physical data `(η,v,x)` from raw data
function mapfro(U)
real(ifft(U[:,1])),real(ifft(U[:,2])),x
end
new( label, f!, mapto, mapfro )
end
A useful tool used in the above was the function Mesh
. It takes as argument a NamedTuple
containing typically a number of points/modes and the (half-)size of the domain and provides the vector of collocations points and Fourier wavenumbers (see this discussion).
The Airy model can now be built as follows
using WaterWaves1D
# include your file
model = Airy((μ=1,L=2π,N=2^8))
build your initial data
The simplest way to build an initial data is to use the function Init
, which takes as argument either
- a function
η
and a functionv
(in this order); - an array of collocation points and two vectors representing
η(x)
andv(x)
(in this order); - a
mesh
(generated withMesh
) and two vectors representingη(mesh.x)
andv(mesh.x)
(in this order).
Some relevant initial data (e.g. travelling waves) are built-in; see the library section. You can build your own in the following lines.
struct Heap <: InitialData
η
v
label :: String
info :: String
function Heap(L)
η=x->exp.(-(L*x).^2)
v=x->zero(x)
init = Init(η,v)
label = "Heap"
info = "Heap of water, with length L=$L."
new( init.η,init.v,label,info )
end
end
The corresponding initial data can now be built as follows
using WaterWaves1D
# include your file
init = Heap(1)
build your time solver
As an example, let us review how the explicit Euler solver, Euler
, is built.
In a dedicated file we write
struct Euler <: TimeSolver
U1 :: Array
label :: String
function Euler( U :: Array; realdata=nothing )
U1 = copy(U)
if realdata==true
U1 = real.(U1);
end
if realdata==false
U1 = complex.(U1);
end
new( U1, "Euler" )
end
end
Here, Euler.U1
is a pre-allocated vector which can be used to speed-up calculations, and Euler-label
is the string "Euler"
, used for future references. The optional keyword argument realdata
allows to specify the type of data which the solver will take as arguments: either complex or real vectors.
We shall now add one method to the function step!
, performing the explicit Euler step: that is replacing a vector U
with U+dt*f(U)
where f
is provided by the model at stake.
export step!
function step!(solver :: Euler,
model :: AbstractModel,
U ,
dt )
solver.U1 .= U # allocate U to U1
model.f!( solver.U1 ) # model.f!(U) replaces its argument U with f(U)
U .+= dt .* solver.U1 # update U
end
access to and manage your data
Once an initial-value problem problem
has been solved (i.e. numerically integrated), the raw data is stored in problem.data.U
, which is an array whose elements correspond (in chronological order) to values at the different computed times, problem.times.ts
.
param = ( c = 1.1, ϵ = 1, μ = 1, N = 2^8, L = 16, T = 1, dt = 0.1)
init = SolitaryWhitham(param)
model = Airy(param)
problem = Problem( model, init, param )
solve!(problem; verbose=false)
problem.data.U
11-element Vector{Matrix{ComplexF64}}:
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -7.002608952054094 + 9.133738484642437e-10im -6.762869498490986 + 8.821034309876048e-10im; … ; 5.494478771973449 + 1.4333285430199028e-9im 5.275526526412349 + 1.376211014282914e-9im; -7.002608952054094 - 9.133737816159736e-10im -6.762869498490985 - 8.821033763093476e-10im]
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -7.001276215168337 + 0.13109975385088155im -6.761582388828596 + 0.1374871831111701im; … ; 5.49044781042491 + 0.197090166103864im 5.27165619666921 + 0.21571490657596468im; -7.001276215168337 - 0.1310997538508815im -6.761582388828595 - 0.13748718311117006im]
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -6.997278511843726 + 0.26214960496561857im -6.757721549808015 + 0.27492203220578515im; … ; 5.478360840660886 + 0.39389114488289im 5.26005088663012 + 0.43111329822862854im; -6.997278511843726 - 0.2621496049656185im -6.757721549808014 - 0.27492203220578515im]
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -6.990617363766434 + 0.3930996714298969im -6.751288451019376 + 0.41225223495202734im; … ; 5.458235597615398 + 0.5901141762054911im 5.240727624500813 + 0.6458791272116842im; -6.990617363766434 - 0.3930996714298968im -6.751288451019375 - 0.41225223495202734im]
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -6.981295306436494 + 0.5239001083980673im -6.74228554115801 + 0.5494255179685713im; … ; 5.4301016105965925 + 0.7854713465158193im 5.2137147628612075 + 0.8596972725473321im; -6.981295306436494 - 0.5239001083980672im -6.742285541158009 - 0.5494255179685713im]
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -6.969315888202691 + 0.6545011279794661im -6.730716247092376 + 0.6863896676040324im; … ; 5.394000159959061 + 0.9796760127160544im 5.179051937064312 + 1.0722540037722788im; -6.969315888202691 - 0.654501127979466im -6.730716247092375 - 0.6863896676040324im]
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -6.954683668911912 + 0.7848530181896827im -6.716584972559662 + 0.8230925498115762im; … ; 5.3499842165341285 + 1.1724432227512511im 5.136790007080184 + 1.2832374412672538im; -6.954683668911912 - 0.7848530181896824im -6.716584972559661 - 0.8230925498115762im]
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -6.93740421817349 + 0.9149061618729335im -6.699897096489541 + 0.9594821299932264im; … ; 5.2981183639069585 + 1.3634901337129568im 5.086990982870261 + 1.4923380138708264im; -6.93740421817349 - 0.9149061618729333im -6.69989709648954 - 0.9594821299932264im]
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -6.917484113239186 + 1.0446110555883417im -6.680658970956735 + 1.0955064928063196im; … ; 5.238478703654544 + 1.5525364268481279im 5.029727933401562 + 1.699248913106074im; -6.917484113239186 - 1.0446110555883414im -6.680658970956734 - 1.0955064928063196im]
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -6.894930936499623 + 1.1739183284529306im -6.658877918763163 + 1.2311138619245658im; … ; 5.171152743683606 + 1.73930471886441im 4.965084879434269 + 1.9036665433536268im; -6.894930936499623 - 1.1739183284529302im -6.658877918763162 - 1.2311138619245658im]
[7.666279752905114 + 0.0im 7.418980406037193 + 0.0im; -6.869753272598125 + 1.302778760934162im -6.634562230650589 + 1.3662526197461948im; … ; 5.096239269832244 + 1.9235209689282842im 4.89315667023999 + 2.105290967310562im; -6.869753272598125 - 1.3027787609341615im -6.634562230650588 - 1.3662526197461948im]
Physical data (say at final computed time) can be reconstructed using the model as follows:
η,v,x = problem.model.mapfro(last(problem.data.U))
([1.3452999698737589e-5, 1.2648564815018193e-5, 1.1943588067414579e-5, 1.1332528165167166e-5, 1.0810582323633744e-5, 1.0373648567388316e-5, 1.0018293546992796e-5, 9.741725616158048e-6, 9.54177295943437e-6, 9.416866602968144e-6 … 2.8821704217644756e-5, 2.6515672474047435e-5, 2.4418100211098977e-5, 2.2512500837691338e-5, 2.0783896157675757e-5, 1.9218698792657085e-5, 1.7804605511088922e-5, 1.653050062875791e-5, 1.538636872894239e-5, 1.436321601783197e-5], [1.3452888212639458e-5, 1.2797168140343818e-5, 1.224207055065818e-5, 1.1783230465495231e-5, 1.1417039625857828e-5, 1.114061811259559e-5, 1.09517916888624e-5, 1.0849074691454685e-5, 1.0831658321136528e-5, 1.0899404252029399e-5 … 2.71342930944201e-5, 2.50326230496669e-5, 2.3127772695689507e-5, 2.140476688537404e-5, 1.985005975398916e-5, 1.845142831639196e-5, 1.7197876421645372e-5, 1.6079548363585883e-5, 1.5087651423169313e-5, 1.4214386778513971e-5], [-16.0, -15.875, -15.75, -15.625, -15.5, -15.375, -15.25, -15.125, -15.0, -14.875 … 14.75, 14.875, 15.0, 15.125, 15.25, 15.375, 15.5, 15.625, 15.75, 15.875])
where η
and v
are respectively the values of the surface deformation and of the derivative of the trace of the velocity potential at the surface, at collocation points x
.
This procedure is carried out by the function solution
, which allows in addition to perform some interpolations (making use of the otherwise helpful function interpolate
).
η,v,x = solution(problem)
([1.3452999698737589e-5, 1.2648564815018193e-5, 1.1943588067414579e-5, 1.1332528165167166e-5, 1.0810582323633744e-5, 1.0373648567388316e-5, 1.0018293546992796e-5, 9.741725616158048e-6, 9.54177295943437e-6, 9.416866602968144e-6 … 2.8821704217644756e-5, 2.6515672474047435e-5, 2.4418100211098977e-5, 2.2512500837691338e-5, 2.0783896157675757e-5, 1.9218698792657085e-5, 1.7804605511088922e-5, 1.653050062875791e-5, 1.538636872894239e-5, 1.436321601783197e-5], [1.3452888212639458e-5, 1.2797168140343818e-5, 1.224207055065818e-5, 1.1783230465495231e-5, 1.1417039625857828e-5, 1.114061811259559e-5, 1.09517916888624e-5, 1.0849074691454685e-5, 1.0831658321136528e-5, 1.0899404252029399e-5 … 2.71342930944201e-5, 2.50326230496669e-5, 2.3127772695689507e-5, 2.140476688537404e-5, 1.985005975398916e-5, 1.845142831639196e-5, 1.7197876421645372e-5, 1.6079548363585883e-5, 1.5087651423169313e-5, 1.4214386778513971e-5], [-16.0, -15.875, -15.75, -15.625, -15.5, -15.375, -15.25, -15.125, -15.0, -14.875 … 14.75, 14.875, 15.0, 15.125, 15.25, 15.375, 15.5, 15.625, 15.75, 15.875], 1.0)
Using (η,v,x)
one can compute other quantities such as the mass, momentum, energy; for instance for the purpose of testing how well these quantities are numerically perserved (when the quantities are first integrals of the considered model). Built-in functions mass
, momentum
, energy
(and massdiff
, momentumdiff
, energydiff
) compute such quantities (or their variation).
plot your data
Once (η,v,x)
is obtained as above, producing plots is as simple as
using Plots
plot(x, [ η v ], label = ["η" "v"])
One can also plot the amplitude of discrete Fourier coefficients (in semi-log scale) as follows:
using FFTW
k=fftshift(Mesh(x).k);fftη=fftshift(fft(η));fftv=fftshift(fft(v));
indices = (fftη .!=0) .& (fftv .!=0 )
plot(k[indices], [ abs.(fftη)[indices] abs.(fftv)[indices] ], yscale=:log10)
Built-in plot recipes provide convenient ways of producing such plots, and animations.
plot(problem;var=[:surface,:velocity,:fourier])
save and load your data
Thanks to the package HDF5.jl
it is possible to save (raw) data to a local file, and then load them for future analyses.
Built-in functions ease the process. Here is a typical example.
Build and solve a problem.
param = ( c = 1.1, ϵ = 1, μ = 1, N = 2^8, L = 16, T = 1, dt = 0.1)
init = SolitaryWhitham(param)
model = Airy(param)
problem = Problem( model, init, param )
solve!(problem; verbose=false)
nothing
[ Info: Converged : relative error 3.472628630637235e-15 in 4 steps
Save the initial data.
x = Mesh(param).x
dump("file_name", x, init)
Save the time-integrated (raw) data.
dump("file_name", problem) # or dump("file_name", problem.data)
Reconstruct the problem, without solving it.
loaded_init = load_init("file_name") # load initial data
new_problem = Problem( model, loaded_init, param ) # re-build problem
load_data!("file_name", new_problem) # incorporate time-integrated data
problem.data == new_problem.data
true
Initial data are functions. The dump
function saves values at some collocation points. Hence the loaded initial data typically differs from the original one by machine epsilon rounding errors.
It is not possible to save models, solvers, or problems. Hence the user needs to store separately the parameters and information required to build them.