Library (index)
WaterWaves1D.AbstractModel
WaterWaves1D.Airy
WaterWaves1D.AkersNicholls
WaterWaves1D.AkersNicholls_fast
WaterWaves1D.Boussinesq
WaterWaves1D.Choi
WaterWaves1D.CnoidalSGN
WaterWaves1D.Data
WaterWaves1D.Euler
WaterWaves1D.EulerSymp
WaterWaves1D.Euler_naive
WaterWaves1D.Init
WaterWaves1D.InitialData
WaterWaves1D.IsobeKakinuma
WaterWaves1D.Matsuno
WaterWaves1D.Matsuno_fast
WaterWaves1D.Mesh
WaterWaves1D.NonHydrostatic
WaterWaves1D.Problem
WaterWaves1D.RK4
WaterWaves1D.RK4_naive
WaterWaves1D.Random
WaterWaves1D.SaintVenant
WaterWaves1D.SerreGreenNaghdi
WaterWaves1D.SolitarySGN
WaterWaves1D.SolitaryWB
WaterWaves1D.SolitaryWGN
WaterWaves1D.SolitaryWhitham
WaterWaves1D.SquareRootDepth
WaterWaves1D.Times
WaterWaves1D.WWn
WaterWaves1D.WaterWaves
WaterWaves1D.Whitham
WaterWaves1D.WhithamBoussinesq
WaterWaves1D.WhithamGreenNaghdi
WaterWaves1D.modifiedMatsuno
WaterWaves1D.relaxedGreenNaghdi
Base.dump
Base.dump
Base.dump
WaterWaves1D.CnoidalWaveSerreGreenNaghdi
WaterWaves1D.SolitaryWaveSerreGreenNaghdi
WaterWaves1D.SolitaryWaveWhitham
WaterWaves1D.SolitaryWaveWhithamBoussinesq
WaterWaves1D.SolitaryWaveWhithamGreenNaghdi
WaterWaves1D.energy
WaterWaves1D.energydiff
WaterWaves1D.interpolate
WaterWaves1D.interpolate
WaterWaves1D.load_data
WaterWaves1D.load_data!
WaterWaves1D.load_init
WaterWaves1D.mass
WaterWaves1D.massdiff
WaterWaves1D.momentum
WaterWaves1D.momentumdiff
WaterWaves1D.random
WaterWaves1D.solution
WaterWaves1D.solve!
WaterWaves1D.solve!
Initial data
WaterWaves1D.CnoidalSGN
— TypeCnoidalSGN(param; P=1)
Build the initial data associated with CnoidalWaveSerreGreenNaghdi(param; P=1)
, of type InitialData
, to be used in initial-value problems Problem(model, initial::InitialData, param)
.
WaterWaves1D.CnoidalWaveSerreGreenNaghdi
— MethodCnoidalWaveSerreGreenNaghdi(param; P=1)
Compute the Serre-Green-Naghdi cnoidal wave with prescribed h₀<h₁<h₂
. h₁
is the minimum, h₂
is the maximum of the wave. As h₀ -> h₁
, the cnoidal wave converges towards the solitary wave. See for instance Gavrilyuk, Nkonga, Shyue and Truskinovsky.
Arguments
param :: NamedTuple
: parameters of the problem containingh₀<h₁<h₂
and dimensionless parametersϵ
andμ
, and number of collocation pointsN
.P :: Int
: (keyword, optional, default = 1) the number of periods of the cnoidal wave in the constructed mesh.
Return values
(η,u,v,mesh,param)
with
η :: Vector{Float64}
: surface deformation;u :: Vector{Float64}
: layer-averaged velocity;v :: Vector{Float64}
: derivative of the trace of the velocity potential at the surface;mesh :: Mesh
: mesh collocation points;param :: NamedTuple
: useful parameters
WaterWaves1D.SolitarySGN
— TypeSolitarySGN(param; x₀=0)
Build the initial data associated with SolitaryWaveSerreGreenNaghdi(param; x₀=0)
, of type InitialData
, to be used in initial-value problems Problem(model, initial::InitialData, param)
.
SolitarySGN(c; ϵ=1,μ=1,x₀=0,N=2^12)
Build the initial data with velocity c
, center x₀
, dimensionless parameters ϵ
and μ
, and number of collocation points N
.
WaterWaves1D.SolitaryWaveSerreGreenNaghdi
— MethodSolitaryWaveSerreGreenNaghdi(param; x₀=0)
Compute the Serre-Green-Naghdi solitary wave with prescribed velocity.
Arguments
param :: NamedTuple
: parameters of the problem containing velocityc
and dimensionless parametersϵ
andμ
, and mesh sizeL
and number of collocation pointsN
;x₀ :: Real
: (keyword, optional, default = 0) center of solitary wave.
Return values
(η,u,v,mesh)
with
η :: Vector{Float64}
: surface deformation;u :: Vector{Float64}
: layer-averaged velocity;v :: Vector{Float64}
: derivative of the trace of the velocity potential at the surface;mesh :: Mesh
: mesh collocation points.
WaterWaves1D.SolitaryWGN
— TypeSolitaryWGN(param; kwargs)
Build the initial data associated with SolitaryWaveWhithamGreenNaghdi(param; kwargs)
, of type InitialData
, to be used in initial-value problems Problem(model, initial::InitialData, param)
.
SolitaryWGN(c; ϵ=1,μ=1,N=2^12,kwargs)
Build the initial data with velocity c
, dimensionless parameters ϵ
and μ
, and number of collocation points N
, and kwargs
the other (optional) keyword arguments as above.
WaterWaves1D.SolitaryWaveWhithamGreenNaghdi
— MethodSolitaryWaveWhithamGreenNaghdi(param; kwargs)
Compute the Whitham-Green-Naghdi solitary wave with prescribed velocity.
Arguments
param :: NamedTuple
: parameters of the problem containing velocityc
and dimensionless parametersϵ
andμ
, and mesh sizeL
and number of collocation pointsN
;
Keywords (optional)
guess :: Vector{Real}
: initial guess for the surface deformation (if not provided, the exact formula for SGN is used);x₀ :: Real
: center of solitary wave (if guess is not provided);SGN :: Bool
: iftrue
computes the Serre-Green-Naghdi (instead of Whitham-Green-Naghdi) solitary wave (considerSolitaryWaveSerreGreenNaghdi
instead);method :: Int
: equation used (between1
and4
);iterative :: Bool
: inverts Jacobian through GMRES iftrue
, LU decomposition iffalse
(default isfalse
);verbose :: Bool
: prints numerical errors at each step iftrue
(default isfalse
);max_iter :: Int
: maximum number of iterations of the Newton algorithm (default is20
);tol :: Real
: relative tolerance measured in ℓ∞ norm (default is1e-10
);ktol :: Real
: tolerance of the Krasny filter (default is0
, i.e. no filtering);gtol :: Real
: relative tolerance of the GMRES algorithm (default is1e-10
);dealias :: Int
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);q :: Real
: Newton algorithm modified withu_{n+1}=q*u_{n+1}+(1-q)*u_n
(default is1
);α :: Real
: addsα
times spectral projection onto the Kernel to the Jacobian (default is0
).
Return values
(η,u,v,mesh)
with
η :: Vector{Float64}
: surface deformation;u :: Vector{Float64}
: layer-averaged velocity;v :: Vector{Float64}
: derivative of the velocity potential at the surface;mesh :: Mesh
: mesh collocation points.
WaterWaves1D.SolitaryWB
— TypeSolitaryWB(param; kwargs)
Build the initial data associated with SolitaryWaveWhithamBoussinesq(param; kwargs)
, of type InitialData
, to be used in initial-value problems Problem(model, initial::InitialData, param)
.
SolitaryWB(c; ϵ=1,μ=1,N=2^12,kwargs)
Build the initial data with velocity c
, dimensionless parameters ϵ
and μ
, and number of collocation points N
, and kwargs
the other (optional) keyword arguments as above.
WaterWaves1D.SolitaryWaveWhithamBoussinesq
— MethodSolitaryWaveWhithamBoussinesq(param; kwargs)
Compute the Whitham-Boussinesq solitary wave with prescribed velocity.
Argument
param :: NamedTuple
: parameters of the problem containing velocityc
and dimensionless parametersϵ
andμ
,
and mesh size L
and number of collocation points N
;
Keywords (optional)
guess :: Vector{Real}
: initial guess for the surface deformation (if not provided, the exact formula for SGN is used);x₀ :: Real
: center of solitary wave (if guess is not provided);α :: Real
: determines the model used (typically1
or1/2
, default is 1);Boussinesq
: iftrue
(default isfalse
), compute the standard Boussinesq system with parametersa
(defaut-1//3
),b=d
(defaut1//3
), andc=0
);iterative :: Bool
: inverts Jacobian through GMRES iftrue
, LU decomposition iffalse
;verbose :: Bool
: prints numerical errors at each step iftrue
;max_iter :: Int
: maximum number of iterations of the Newton algorithm;tol :: Real
: general tolerance (default is1e-10
);ktol :: Real
: tolerance of the Krasny filter (default is0
, i.e. no filtering);gtol :: Real
: relative tolerance of the GMRES algorithm (default is1e-10
);dealias :: Int
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);q :: Real
: Newton algorithm modified with
u_{n+1}=q*(u_n+du)+(1-q)*u_n
(default is 1
);
β :: Real
: addsβ
times spectral projection onto the Kernel to the Jacobian.β
Return values
(η,v,mesh)
with
η :: Vector{Float64}
: surface deformation;v :: Vector{Float64}
: velocity (derivative of the trace of the velocity potential at the surface);mesh :: Mesh
: mesh collocation points.
WaterWaves1D.SolitaryWhitham
— TypeSolitaryWhitham(param; kwargs)
Build the initial data associated with SolitaryWaveWhitham(param; kwargs)
, of type InitialData
, to be used in initial-value problems Problem(model, initial::InitialData, param)
.
SolitaryWhitham(c; ϵ=1,μ=1,N=2^12,kwargs)
Build the initial data with velocity c
, dimensionless parameters ϵ
and μ
, and number of collocation points N
, and kwargs
the other (optional) keyword arguments as above.
WaterWaves1D.SolitaryWaveWhitham
— MethodSolitaryWaveWhitham(param; kwargs)
Compute the Whitham solitary wave with prescribed velocity.
Argument
param :: NamedTuple
: parameters of the problem containing velocityc
and dimensionless parametersϵ
andμ
, and mesh sizeL
and number of collocation pointsN
;
Keywords (optional)
guess :: Vector{Real}
: initial guess for the surface deformation (if not provided, the exact formula for KdV is used);x₀ :: Real
: center of solitary wave (if guess is not provided);iterative :: Bool
: inverts Jacobian through GMRES iftrue
, LU decomposition iffalse
;verbose :: Bool
: prints numerical errors at each step iftrue
;max_iter :: Int
: maximum number of iterations of the Newton algorithm;tol :: Real
: general tolerance (default is1e-10
);ktol :: Real
: tolerance of the Krasny filter (default is0
, i.e. no filtering);gtol :: Real
: relative tolerance of the GMRES algorithm (default is1e-10
);dealias :: Int
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);q :: Real
: Newton algorithm modified with
u_{n+1}=q*u_{n+1}+(1-q)*u_n
(default is 1
);
α :: Real
: addsα
times spectral projection onto the Kernel to the Jacobian;KdV :: Bool
: iftrue
computes the KdV (instead of Whitham) solitary wave.
Return values
(u,mesh)
with
u :: Vector{Float64}
: the solution;mesh :: Mesh
: mesh collocation points.
WaterWaves1D.Random
— TypeRandom(param;args)
Randomly generated initial data, based on provided (optional arguments) :
L
is the typical wavelength (default isL=1
),s
is the (real) Sobolev index regularity (default iss=∞
),λ
is the length of spatial localization (default isλ=∞
, no localization),a
is the couple of amplitudes of the surface deformation, and velocity (default isa=(1,1)
).
Return an initial data init::InitialData
, to be used in initial-value problems Problem(model, init, param)
, with init.η
and init.v
two samples of random(x;L,s,λ)
, where x
is the set of collocation points generated by Mesh(param)
.
WaterWaves1D.random
— Methodrandom(x;args)
Randomly generate a vector of values at collocation points x
, based on provided (optional arguments) :
L
is the typical wavelength (default isL=1
),s
is the (real) Sobolev index regularity (default iss=∞
, smooth data),λ
is the length of spatial localization (default isλ=∞
, no localization),a
is the amplitude of the returned vector (default isa=1
).
The vector is generated through randomly chosen Fourier coefficients, multiplied with weigth w=10^(-|k|L/(2π))
if s=∞
, or w=1/(1+9(|k|L/(2π))^(s+1/2))
otherwise. If λ≠∞
, the function in spatial variables is multiplied by exp(-|x/λ|^2)
, and in any case normalized to have maximum absolute value 1.
Models
WaterWaves1D.AbstractModel
— TypeAbstract type whose subtypes are the models from which initial-value problems can be built, through Problem( model :: AbstractModel, initial :: InitialData, param :: NamedTuple )
WaterWaves1D.Airy
— TypeAiry(param;mesh,label)
Define an object of type AbstractModel
in view of solving the initial-value problem for the linear (Airy) water waves equations.
Arguments
param::NamedTuple
must contain- the shallowness parameter
μ
(set the infinite-layer case ifμ=Inf
); - optionally,
ν
the shallow/deep water scaling factor. By default,ν=1
ifμ≦1
andν=1/√μ
otherwise; - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided.
- the shallowness parameter
mesh :: Mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
.label :: String
: a label for future references (default is"linear (Airy)"
).
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
Airy.f!
to be called in explicit time-integration solvers; - a function
Airy.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
Airy.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
.
WaterWaves1D.AkersNicholls
— TypeAkersNicholls(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the quadratic deep-water model proposed by Akers and Nicholls and Cheng, Granero-Belinchón, Shkoller and Milewski
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - optionally,
ν
the shallow/deep water scaling factor. By default,ν=1
ifμ≦1
andν=1/√μ
otherwise. Set the infinite-layer case ifν=0
, orμ=Inf
. - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;IL
: Set the infinite-layer case ifIL=true
(orμ=Inf
, orν=0
), in which caseϵ
is the steepness parameter. Default isfalse
.dealias
: dealiasing with1/3
Orlicz rule iftrue
or no dealiasing iffalse
(by default);label
: a label for future references (default is"deep quadratic"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
AkersNicholls.f!
to be called in explicit time-integration solvers; - a function
AkersNicholls.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
AkersNicholls.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
.
WaterWaves1D.AkersNicholls_fast
— TypeAkersNicholls_fast(param;kwargs)
Same as AkersNicholls
, but faster.
WaterWaves1D.Boussinesq
— TypeBoussinesq(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for abcd
-Boussinesq models (with b=d
and c=0
). See Bona, Chen, and Saut
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
- two parameters
a
(default is-1/3
) andb
(default is+1/3
) which determine the model solved. You needa+2*b=1/3
for validity as a long wave model (without surface tension). mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;ktol
: tolerance of the low-pass Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"Boussinesq"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
Boussinesq.f!
to be called in explicit time-integration solvers; - a function
Boussinesq.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
Boussinesq.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
.
WaterWaves1D.Choi
— TypeChoi(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the model proposed by Choi.
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
M∈{0,1,2}
: the order of the model.M=2
by default;reg
: applies a regularization operator.reg=false
by default.mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;iterative
: solve the elliptic problem through GMRES iftrue
, LU decomposition iffalse
(default istrue
);precond
: use a (left) preconditioner for GMRES iftrue
(default), chooseprecond
as the preconditioner if provided;gtol
: relative tolerance of the GMRES algorithm (default is1e-14
);restart
: the corresponding option of the GMRES algorithm (default is100
);maxiter
: the corresponding option of GMRES (default isnothing
);ktol
: tolerance of the Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"Choi-N"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
Choi.f!
to be called in explicit time-integration solvers; - a function
Choi.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
Choi.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, where
- `η` is the values of surface deformation at collocation points `x`;
- `v` is the derivative of the trace of the velocity potential;
- additionally, a handy function
Choi.mapfrofull
which from data matrix returns the Tuple of real vectors(η,v,u,vᵦ)
, where
- `v` is the derivative of the trace of the velocity potential;
- `u` corresponds to the layer-averaged horizontal velocity;
- `vᵦ` corresponds to the horizontal velocity at the bottom.
WaterWaves1D.IsobeKakinuma
— TypeIsobeKakinuma(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the Isobe-Kakinuma model proposed by Isobe.
Argument
param
is of type NamedTuple
(or a collection NamedTuple
s) of and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;iterative
: solve the elliptic problem through GMRES iftrue
, LU decomposition iffalse
(default istrue
);precond
: use a (left) preconditioner for GMRES iftrue
(default), chooseprecond
as the preconditioner if provided;gtol
: relative tolerance of the GMRES algorithm (default is1e-14
);restart
: the corresponding option of the GMRES algorithm (default is100
);maxiter
: the corresponding option of GMRES (default isnothing
);ktol
: tolerance of the Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"Isobe-Kakinuma"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
IsobeKakinuma.f!
to be called in explicit time-integration solvers; - a function
IsobeKakinuma.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
IsobeKakinuma.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
.
WaterWaves1D.Matsuno
— TypeMatsuno(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the quadratic deep-water model proposed by Matsuno.
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - optionally,
ν
the shallow/deep water scaling factor. By default,ν=1
ifμ≦1
andν=1/√μ
otherwise. Set the infinite-layer case ifν=0
, orμ=Inf
. - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
IL
: Set the infinite-layer case ifIL=true
(orμ=Inf
, orν=0
), in which caseϵ
is the steepness parameter. Default isfalse
.mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;dealias
: dealiasing with1/3
Orlicz rule iftrue
or no dealiasing iffalse
(by default);label
: a label for future references (default is"Matsuno"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
Matsuno.f!
to be called in explicit time-integration solvers; - a function
Matsuno.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
Matsuno.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
.
WaterWaves1D.Matsuno_fast
— TypeMatsuno_fast(param;kwargs)
Same as Matsuno
, but faster.
WaterWaves1D.NonHydrostatic
— TypeNonHydrostatic(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the "Non-hydrostatic" model proposed by Bristeau, Mangeney, Sainte-Marie and Seguin
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;iterative
: solve the elliptic problem through GMRES iftrue
, LU decomposition iffalse
(default istrue
);precond
: use a (left) preconditioner for GMRES iftrue
(default), chooseprecond
as the preconditioner if provided;gtol
: relative tolerance of the GMRES algorithm (default is1e-14
);restart
: the corresponding option of the GMRES algorithm (default is100
);maxiter
: the corresponding option of GMRES (default isnothing
);ktol
: tolerance of the Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"non-hydrostatic"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
NonHydrostatic.f!
to be called in explicit time-integration solvers; - a function
NonHydrostatic.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
NonHydrostatic.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
;
- additionally, a handy function
NonHydrostatic.mapfrofull
which from data matrix returns the Tuple of real vectors(η,v,u)
, whereu
corresponds to the layer-averaged velocity.
WaterWaves1D.SaintVenant
— TypeSaintVenant(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for Saint-Venant (or shallow water) model.
Argument
param
is of type NamedTuple
and must contain
- the dimensionless parameter
ϵ
(nonlinearity); - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;ktol
: tolerance of the low-pass Krasny filter (default is0
, i.e. no filtering);dealias
: no dealisasing if set to0
orfalse
(default), standard 3/2 Orlicz rule if set to1
ortrue
, otherwise the value sets additionnally a maximal slope of the dealiasing symbol (2/dealias
models are affected);label
: a label for future references (default is"Saint-Venant"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
SaintVenant.f!
to be called in explicit time-integration solvers; - a function
SaintVenant.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
SaintVenant.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
.
WaterWaves1D.SerreGreenNaghdi
— TypeSerreGreenNaghdi(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the Serre-Green-Naghdi model (Serre, Su and Gardner, Green and Naghdi).
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;iterative
: solve the elliptic problem through GMRES iftrue
, LU decomposition iffalse
(default istrue
);precond
: use a (left) preconditioner for GMRES iftrue
(default), chooseprecond
as the preconditioner if provided;gtol
: relative tolerance of the GMRES algorithm (default is1e-14
);restart
: the corresponding option of the GMRES algorithm (default is100
);maxiter
: the corresponding option of GMRES (default isnothing
);ktol
: tolerance of the Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"Green-Naghdi"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
SerreGreenNaghdi.f!
to be called in explicit time-integration solvers; - a function
SerreGreenNaghdi.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
SerreGreenNaghdi.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
;
- additionally, a handy function
SerreGreenNaghdi.mapfrofull
which from data matrix returns the Tuple of real vectors(η,v,u)
, whereu
corresponds to the layer-averaged velocity.
WaterWaves1D.SquareRootDepth
— TypeSquareRootDepth(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the "√D" model proposed by Cotter, Holm and Percival
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;iterative
: solve the elliptic problem through GMRES iftrue
, LU decomposition iffalse
(default istrue
);precond
: use a (left) preconditioner for GMRES iftrue
(default), chooseprecond
as the preconditioner if provided;gtol
: relative tolerance of the GMRES algorithm (default is1e-14
);restart
: the corresponding option of the GMRES algorithm (default is100
);maxiter
: the corresponding option of GMRES (default isnothing
);ktol
: tolerance of the Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"square-root depth"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
SquareRootDepth.f!
to be called in explicit time-integration solvers; - a function
SquareRootDepth.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
SquareRootDepth.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
;
- additionally, a handy function
SquareRootDepth.mapfrofull
which from data matrix returns the Tuple of real vectors(η,v,u)
, whereu
corresponds to the layer-averaged velocity.
WaterWaves1D.WWn
— TypeWWn(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the water waves expansion proposed by Dommermuth and Yue, West et al., Craig and Sulem (see also the account by Choi) with the "rectification" method proposed by Duchêne and Melinand.
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - optionally,
ν
the shallow/deep water scaling factor. By default,ν=1
ifμ≦1
andν=1/√μ
otherwise. Set the infinite-layer case ifν=0
, orμ=Inf
. - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
IL
: Set the infinite-layer case ifIL=true
(orμ=Inf
, orν=0
), in which caseϵ
is the steepness parameter. Default isfalse
.n :: Int
: the order of the expansion; linear system if1
, quadratic if2
, cubic if3
, quartic if4
(default and other values yield2
);δ
andm
: parameters of the rectifier operator, set ask->min(1,|δ*k|^m)
ork->min(1,|δ*k|^m[1]*exp(1-|δ*k|^m[2]))
ifm
is a couple
(by default is δ=0
, i.e. no regularization and m=-1
. Notice m=-Inf
and δ>0
yields a cut-off filter);
mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;ktol
: tolerance of the low-pass Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"WWn"
withn
the order of the expansion);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
WWn.f!
to be called in explicit time-integration solvers (alsoWWn.f1!
andWWn.f2!
for the symplectic Euler solver); - a function
WWn.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
WWn.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, where
- `η` is the values of surface deformation at collocation points `x`;
- `v` is the derivative of the trace of the velocity potential at `x`.
WaterWaves1D.WaterWaves
— TypeWaterWaves(param; kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the water waves system (via conformal mapping, see Zakharov, Dyachenko and Vasilyev).
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - optionally,
ν
the shallow/deep water scaling factor. By default,ν=1
ifμ≦1
andν=1/√μ
otherwise. Set the infinite-layer case ifν=0
, orμ=Inf
. - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
IL
: Set the infinite-layer case ifIL=true
, in which caseϵ
is the steepness parameter. Default isfalse
.mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;method ∈ {1,2,3}
: method used to initialize the conformal mapping, as a fix-point problemF(u)=u
- if
method == 1
, use standard contraction fix-point iteration; - if
method == 2
, use Newton algorithm with GMRES iterative solver to invert the Jacobian; - if
method == 3
, use Newton algorithm with direct solver to invert the Jacobian;
- if
tol
: (relative) tolerance of the fix-point algorithm (default is1e-16
);maxiter
: the maximal number of iteration in the fix-point algorithm (default is100
);ktol
: tolerance of the low-pass Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"water waves"
);verbose
: prints information iftrue
(default istrue
).
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
WaterWaves.f!
to be called in the explicit time-integration solver (alsoWaterWaves.f1!
andWaterWaves.f2!
for the symplectic Euler solver); - a function
WaterWaves.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
WaterWaves.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, wherex
is a vector of collocation points (non-regularly spaced);
- `η` is the values of surface deformation at collocation points `x`;
- `v` is the derivative of the trace of the velocity potential at `x`.
WaterWaves1D.Whitham
— TypeWhitham(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for two uncoupled Whitham equations, following Emerald.
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
KdV
: iftrue
(default isfalse
), compute the standard KdV equations instead (seeKdV(param;kwargs)
);mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;ktol
: tolerance of the low-pass Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"Whitham-Boussinesq"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
Whitham.f!
to be called in explicit time-integration solvers; - a function
Whitham.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed. - a function
Whitham.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, where
- `η` is the values of surface deformation at collocation points `x`;
- `v` is the derivative of the trace of the velocity potential at `x`.
WaterWaves1D.WhithamBoussinesq
— TypeWhithamBoussinesq(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for a Boussinesq-type model with full-dispersion property.
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
Boussinesq
: iftrue
(default isfalse
), compute the standard Boussinesq system instead (seeBoussinesq(param;kwargs)
);- a parameter
α
which determines the model solved:- If
α = 1
(default), then the model has been introduced in Dinvay, Dutykh and Kalisch; - If
α = 1/2
, then the model is a quasilinear version; - If
α < 1/2
, then expect instabilities stemming from ill-posedness of the model.
- If
mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;ktol
: tolerance of the low-pass Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"Whitham-Boussinesq"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
WhithamBoussinesq.f!
to be called in explicit time-integration solvers; - a function
WhithamBoussinesq.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed. - a function
WhithamBoussinesq.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, where
- `η` is the values of surface deformation at collocation points `x`;
- `v` is the derivative of the trace of the velocity potential at `x`.
WaterWaves1D.WhithamGreenNaghdi
— TypeWhithamGreenNaghdi(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the fully dispersive Green-Naghdi model proposed by Duchêne, Israwi and Talhouk.
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
SGN
: iftrue
(default isfalse
), compute the Serre-Green-Naghdi (SGN) instead of Whitham-Green-Naghdi (WGN) system (seeSerreGreenNaghdi(param;kwargs)
);mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;iterative
: solve the elliptic problem through GMRES iftrue
, LU decomposition iffalse
(default istrue
);precond
: use a (left) preconditioner for GMRES iftrue
(default), chooseprecond
as the preconditioner if provided;gtol
: relative tolerance of the GMRES algorithm (default is1e-14
);restart
: the corresponding option of the GMRES algorithm (default is100
);maxiter
: the corresponding option of GMRES (default isnothing
);ktol
: tolerance of the Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"Whitham-Green-Naghdi"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
WhithamGreenNaghdi.f!
to be called in explicit time-integration solvers; - a function
WhithamGreenNaghdi.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
WhithamGreenNaghdi.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, where
- `η` is the values of surface deformation at collocation points `x`;
- `v` is the derivative of the trace of the velocity potential at `x`;
- additionally, a handy function
WhithamGreenNaghdi.mapfrofull
which from data matrix returns the Tuple of real vectors(η,v,u)
, whereu
corresponds to the layer-averaged velocity.
WaterWaves1D.modifiedMatsuno
— TypemodifiedMatsuno(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for the modified Matsuno model proposed by Duchêne and Melinand.
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - optionally,
ν
the shallow/deep water scaling factor. By default,ν=1
ifμ≦1
andν=1/√μ
otherwise. Set the infinite-layer case ifν=0
, orμ=Inf
. - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
IL
: Set the infinite-layer case ifIL=true
(orμ=Inf
, orν=0
), in which caseϵ
is the steepness parameter. Default isfalse
.mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;ktol
: tolerance of the low-pass Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"modified Matsuno"
);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
modifiedMatsuno.f!
to be called in explicit time-integration solvers; - a function
modifiedMatsuno.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
modifiedMatsuno.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, whereη
is the values of surface deformation at collocation pointsx
;v
is the derivative of the trace of the velocity potential atx
.
WaterWaves1D.relaxedGreenNaghdi
— TyperelaxedGreenNaghdi(param;kwargs)
Define an object of type AbstractModel
in view of solving the initial-value problem for a relaxed Green-Naghdi model proposed by N. Favrie and S. Gavrilyuk or C. Escalante, M. Dumbser and M. Castro and G. Richard.
Argument
param
is of type NamedTuple
and must contain
- dimensionless parameters
ϵ
(nonlinearity) andμ
(dispersion); - the relaxation parameter
a
; - numerical parameters to construct the mesh of collocation points, if
mesh
is not provided as a keyword argument.
Optional keyword arguments
FG
: iftrue
(default isfalse
), compute the Favrie-Gavrilyuk model, otherwise compute the Escalante-Dumbser-Castro model;id
:∈{0,1,2}
and represent the level of preparation of the initial data (default is1
);mesh
: the mesh of collocation points. By default,mesh = Mesh(param)
;iterative
: solve the elliptic problem (to construct initial data) through GMRES iftrue
, LU decomposition iffalse
(default istrue
);precond
: use a (left) preconditioner for GMRES iftrue
(default), chooseprecond
as the preconditioner if provided;gtol
: relative tolerance of the GMRES algorithm (default is1e-14
);restart
: the corresponding option of the GMRES algorithm (default is100
);maxiter
: the corresponding option of GMRES (default isnothing
);ktol
: tolerance of the Krasny filter (default is0
, i.e. no filtering);dealias
: dealiasing with Orlicz rule1-dealias/(dealias+2)
(default is0
, i.e. no dealiasing);label
: a label for future references (default is"Favrie-Gavrilyuk"
ifFG==true
,"Escalante-Dumbser-Castro"
otherwise);
Return values
Generate necessary ingredients for solving an initial-value problem via solve!
:
- a function
relaxedGreenNaghdi.f!
to be called in explicit time-integration solvers; - a function
relaxedGreenNaghdi.mapto
which from(η,v)
of typeInitialData
provides the raw data matrix on which computations are to be executed; - a function
relaxedGreenNaghdi.mapfro
which from such data matrix returns the Tuple of real vectors(η,v,x)
, where
- `η` is the values of surface deformation at collocation points `x`;
- `v` is the derivative of the trace of the velocity potential at `x`;
- additionally, a handy function
relaxedGreenNaghdi.mapfrofull
which from data matrix returns the Tuple of real vectors(η,v,u,p,w)
, where
- `u` corresponds to the layer-averaged horizontal velocity.
- `p` corresponds to the relaxed (artificial) layer-averaged non-hydrostatic pressure;
- `w` corresponds to the relaxed (artificial) layer-averaged vertical velocity.
Solvers
WaterWaves1D.Euler
— TypeEuler(arguments;realdata)
Explicit Euler solver.
Construct an object of type TimeSolver
to be used in Problem(model, initial, param; solver::TimeSolver)
Arguments can be either
- an object of type
AbstractModel
; - an
Array
of size(N,datasize)
whereN
is the number of collocation points anddatasize
the number of equations solved; (param,datasize)
whereparam is a
NamedTuplecontaining a key
N, and
datasizea integer (optional, by default
datasize=2`).
The keyword argument realdata
is optional, and determines whether pre-allocated vectors are real- or complex-valued. By default, they are either determined by the model or the type of the array in case 0.
and 1.
, complex-valued in case 2.
.
WaterWaves1D.EulerSymp
— TypeEulerSymp(arguments;Niter,implicit,realdata)
Symplectic Euler solver. The implicit Euler method is first used on one equation, then the explicit Euler method is used on the second one. The implicit equation is solved via Neumann iteration
Construct an object of type TimeSolver
to be used in Problem(model, initial, param; solver::TimeSolver)
Arguments can be either
- an object of type
AbstractModel
; - an
Array
of size(N,2)
whereN
is the number of collocation points; - a
NamedTuple
containing a keyN
.
The keyword argument Niter
(optional, defaut value = 10) determines the number of steps in the Neumann iteration solver of the implicit step. The keyword argument implicit
(optional, defaut value = 1) determines which equation is implicit (must be 1
or 2
). The keyword argument realdata
is optional, and determines whether pre-allocated vectors are real- or complex-valued. By default, they are either determined by the model or the type of the array in case 0.
and 1.
, complex-valued in case 2.
.
WaterWaves1D.Euler_naive
— TypeEuler_naive()
Runge-Kutta fourth order solver.
A naive version of Euler
, without argument since no pre-allocation is performed.
WaterWaves1D.RK4
— TypeRK4(arguments;realdata)
Explicit Runge-Kutta fourth order solver.
Construct an object of type TimeSolver
to be used in Problem(model, initial, param; solver::TimeSolver)
Arguments can be either
- an object of type
AbstractModel
, typically the model you will solve with the solver; - an
Array
of size(N,datasize)
whereN
is the number of collocation points anddatasize
the number of equations solved; (param,datasize)
whereparam
is aNamedTuple
containing a keyN
, anddatasize
a integer (optional, by defaultdatasize=2
).
The keyword argument realdata
is optional, and determines whether pre-allocated vectors are real- or complex-valued. By default, they are either determined by the model or the type of the array in case 0.
and 1.
, complex-valued in case 2.
.
WaterWaves1D.RK4_naive
— TypeRK4_naive()
Runge-Kutta fourth order solver.
A naive version of RK4
, without argument since no pre-allocation is performed.
Structures
WaterWaves1D.Init
— TypeInit(data ; fast, label)
Generate an initial data to be used in the function Problem
.
data
should contain either
- a function
η
and a functionv
(in this order); - a
mesh
(of typeMesh
) and two vectors representingη(mesh.x)
andv(mesh.x)
(in this order); - an array of collocation points and two vectors representing
η(x)
andv(x)
(in this order).
In the last two cases, an optional keyword argument fast
can be set to true
, (default is false
), in which case the algorithm is faster and uses less allocations, but is less precise.
In the last case, the collocation points must be regularly spaced, otherwise an ErrorException
is raised.
If the keyword label::String
(used to display information to the output stream) is not provided, then it is set to the "user-defined"
.
WaterWaves1D.InitialData
— TypeAbstract type defining initial data from which initial-value problems can be built, through Problem( model :: AbstractModel, initial :: InitialData, param :: NamedTuple )
WaterWaves1D.Mesh
— TypeMesh(args)
Construct a periodic mesh of N
collocation points regularly spaced between xmin
(included) and xmax
(excluded), and associated Fourier modes.
Arguments
Can be either
- a
NamedTuple
containingN
andxmin
andxmax
; or - a
NamedTuple
containingN
andL
(in which casexmin=-L
andxmax=L
); or - a vector of regularly spaced collocation points.
Return values
m=Mesh(args)
is of parametric type and offers with
m.N
: number of collocation points and Fourier modes;m.xmin
: minimum of the mesh (included in the vector of collocation points);m.xmax
: maximum of the mesh (excluded in the vector of collocation points);m.dx
: distance between two collocation points;m.x
: the vector of collocation points;m.k
: the vector of wavenumbers;m.kmin
: minimum of wavenumbers (included in the vector of wavenumbers);m.kmax
: maximum of wavenumbers (included in the vector of wavenumbers);m.dk
: distance between two Fourier modes.
WaterWaves1D.Problem
— TypeProblem( model, initial, param ; solver, label)
Build an initial-value problem which can then be solved (i.e. integrated in time) through solve!( problem )
Arguments
model :: AbstractModel
, the system of equation solved.
May be built, e.g., by WaterWaves(param)
;
initial :: InitialData
, the initial data.
May be buit, e.g., by Init(η,v)
where η
is the surface deformation and v
the derivative of the trace of the velocity potential at the surface;
times :: Times
is the time grid, and may be built using the functionTimes
.
Alternatively one can simply provide a NamedTuple
with - T
, the final time of integration - dt
, the timestep - optionally, Ns
the number of computed data or ns
for storing data every ns
computation steps (by default, every computed data is stored).
Optional keyword arguments
solver :: TimeSolver
, the solver for time integration (default is explicit Runge-Kutta fourth order solver).
May be built, e.g., by RK4(model)
or RK4_naive()
.
label :: String
is used in future references (e.g.plot_solution
).- Information are not printed if
verbose = false
(default istrue
).
WaterWaves1D.solve!
— Methodsolve!( problems; verbose=true )
Solve (i.e. integrate in time) a collection of initial-value problems.
The argument problems
should be a collection (list, array...) of elements of type Problem
.
Information are not printed if keyword argument verbose = false
(default is true
).
WaterWaves1D.solve!
— Methodsolve!( problem :: Problem; verbose=true )
Solve (i.e. integrate in time) an initial-value problem
The argument problem
should be of type Problem
. It may be buit, e.g., by Problem(model, initial, param)
Information are not printed if keyword argument verbose = false
(default is true
).
WaterWaves1D.Data
— TypeData( m :: Matrix )
Data structure to store the solution of an initial-value problem along time.
data=Data(m)
is of parametric type and offers
data.U
, a 1-element vector with a copy of the matrixm
;(data.datalength,data.datasize)=size(m)
wheredatalength
is the number of computed modes, anddatasize
the number of involved equations, typically 2.
WaterWaves1D.Times
— TypeTimes(param; ns, Ns)
Constructs a mesh of times, to be used in initial-value problems (see Problem
).
Arguments
param
is either
- a
NamedTuple
containingdt
the timestep andT
the final time of comuptation; or - a vector of computed times.
Optional keyword arguments
ns
: data are stored everyns
computations (optional, default = 1).Ns
:Ns
data (in addition to the initial datum) are stored (optional, by default `floor( tfin/dt)).
If both Ns
and ns
are given, Ns
overrules ns
.
Return values
t=Times(args)
is of parametric type and offers
t.Nc
: number of computed times (including initial datum);t.Ns
: number of stored times (including initial datum);t.ns
: number of computed times between two stored times;t.tfin
: the final time;t.dt
: the timestep;t.tc
: the vector of computed times;t.ts
: the vector of stored times.
Tools
WaterWaves1D.energy
— Methodenergy(pb::Problem;T)
Compute the excess of mass of a solved initial-value problem pb
at a given time T
.
Keyword argument T
is optional, the last computed time is used by default.
WaterWaves1D.energydiff
— Methodenergydiff(pb::Problem;T,rel)
Compute the difference of energy of a solved initial-value problem pb
between given time T
and initial time.
Keyword argument T
is optional, the last computed time is used by default.
If keyword argument rel=true
(default is false), then compute the relative difference (with initial value as reference).
WaterWaves1D.interpolate
— Methodinterpolate(mesh,vector,x;fast)
Interpolate a vector vector
of values on a uniform collocation grid defined by mesh
, on collocation points given by x
.
If the collocation points x
are regularly spaced and the optional keyword argument fast
is set to true
(default is false
), then the algorithm is faster and uses less allocations, but is less precise.
Return the vector of values on collocation points.
WaterWaves1D.interpolate
— Methodinterpolate(mesh,vector;n=2^3)
Interpolate a vector vector
of values on a uniform collocation grid defined by mesh
.
Return (new_mesh,new_vector)
a new uniform mesh with n
times as many values, and the vector of values at these collocation points.
WaterWaves1D.mass
— Methodmass(pb::Problem;T)
Compute the excess of mass of a solved initial-value problem pb
at a given time T
.
Keyword argument T
is optional, the last computed time is used by default.
WaterWaves1D.massdiff
— Methodmassdiff(pb::Problem;T,rel)
Compute the difference of excess of mass of a solved initial-value problem pb
between given time T
and initial time.
Keyword argument T
is optional, the last computed time is used by default.
If keyword argument rel=true
(default is false), then compute the relative difference (with initial value as reference).
WaterWaves1D.momentum
— Methodmomentum(pb::Problem;T)
Compute the horizontal impulse of a solved initial-value problem pb
at a given time T
.
Keyword argument T
is optional, the last computed time is used by default.
WaterWaves1D.momentumdiff
— Methodmomentumdiff(pb::Problem;T,rel)
Compute the difference of horizontal impulse of a solved initial-value problem pb
between given time T
and initial time.
Keyword argument T
is optional, the last computed time is used by default.
If keyword argument rel=true
(default is false), then compute the relative difference (with initial value as reference).
WaterWaves1D.solution
— Methodsolution(pb::Problem;T,x,interpolation)
Give the solution of a solved initial-value problem at a given time T
.
Arguments
- Argument
pb
is of typeProblem
. - Keyword argument
T
is optional, the last computed time is returned by default. - Keyword argument
x
is optional, if provided the solution is interpolated to the collocation vectorx
. - Keyword argument
interpolation
is optional, if an integer is provided the solution is interpolated on as many collocation points (iftrue
, then the default value2^3
is chosen). - Keyword argument
raw
is optional, if set totrue
then(U,t)
withU
the raw data andt
the time is returned (default isfalse
).
Return values
Return (η,v,x,t)
where
η
is the surface deformation at collocation points;v
is the tangential velocity (derivative of the trace of the velocity potential) at collocation points;x
is the vector of collocation points;t
the time (first computed time greater or equal to providedT
).
Load and save
Base.dump
— Methoddump(file_name :: String, data :: Data)
Save data
to the file with name file_name
(and extension ".h5").
Base.dump
— Methoddump(file_name :: String, problem :: Problem)
Save the data of problem
to the file with name file_name
(and extension ".h5").
Base.dump
— Methoddump(file_name :: String, x::Vector, init::InitialData )
Save the values of the initial data init
at collocation points x
to the file with name file_name
(and extension ".h5").
WaterWaves1D.load_data!
— Methodload_data!(file_name :: String, problem :: Problem)
Fills problem
with raw data extracted from the file with name file_name
(and extension ".h5").
WaterWaves1D.load_data
— Methodload_data(file_name :: String)
Load data from the file with name file_name
(and extension ".h5").
WaterWaves1D.load_init
— Methodload_init(file_name :: String; fast = false)
Load initial data from the file with name file_name
(and extension ".h5").
Keyword argument fast
is optional (default is false
), and corresponds the the keyword argument of the function Init
.
Return an object of type InitialData
.