Library (index)

Initial data

    WaterWaves1D.CnoidalSGNType
    CnoidalSGN(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).

    source
    WaterWaves1D.CnoidalWaveSerreGreenNaghdiMethod
    CnoidalWaveSerreGreenNaghdi(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 containing h₀<h₁<h₂ and dimensionless parameters ϵ and μ, and number of collocation points N.
    • 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
    source
    WaterWaves1D.SolitarySGNType
    SolitarySGN(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.

    source
    WaterWaves1D.SolitaryWaveSerreGreenNaghdiMethod
    SolitaryWaveSerreGreenNaghdi(param; x₀=0)

    Compute the Serre-Green-Naghdi solitary wave with prescribed velocity.

    Arguments

    • param :: NamedTuple: parameters of the problem containing velocity c and dimensionless parameters ϵ and μ, and mesh size L and number of collocation points N;
    • 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.
    source
    WaterWaves1D.SolitaryWGNType
    SolitaryWGN(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.

    source
    WaterWaves1D.SolitaryWaveWhithamGreenNaghdiMethod
    SolitaryWaveWhithamGreenNaghdi(param; kwargs)

    Compute the Whitham-Green-Naghdi solitary wave with prescribed velocity.

    Arguments

    • param :: NamedTuple: parameters of the problem containing velocity c 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);
    • SGN :: Bool: if true computes the Serre-Green-Naghdi (instead of Whitham-Green-Naghdi) solitary wave (consider SolitaryWaveSerreGreenNaghdi instead);
    • method :: Int: equation used (between 1 and 4);
    • iterative :: Bool: inverts Jacobian through GMRES if true, LU decomposition if false (default is false);
    • verbose :: Bool: prints numerical errors at each step if true (default is false);
    • max_iter :: Int: maximum number of iterations of the Newton algorithm (default is 20);
    • tol :: Real: relative tolerance measured in ℓ∞ norm (default is 1e-10);
    • ktol :: Real: tolerance of the Krasny filter (default is 0, i.e. no filtering);
    • gtol :: Real: relative tolerance of the GMRES algorithm (default is 1e-10);
    • dealias :: Int: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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 (default is 0).

    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.
    source
    WaterWaves1D.SolitaryWBType
    SolitaryWB(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.

    source
    WaterWaves1D.SolitaryWaveWhithamBoussinesqMethod
    SolitaryWaveWhithamBoussinesq(param; kwargs)

    Compute the Whitham-Boussinesq solitary wave with prescribed velocity.

    Argument

    • param :: NamedTuple: parameters of the problem containing velocity c 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 (typically 1 or 1/2, default is 1);
    • Boussinesq: if true (default is false), compute the standard Boussinesq system with parameters a (defaut -1//3), b=d (defaut 1//3), and c=0);
    • iterative :: Bool: inverts Jacobian through GMRES if true, LU decomposition if false;
    • verbose :: Bool: prints numerical errors at each step if true;
    • max_iter :: Int: maximum number of iterations of the Newton algorithm;
    • tol :: Real: general tolerance (default is 1e-10);
    • ktol :: Real: tolerance of the Krasny filter (default is 0, i.e. no filtering);
    • gtol :: Real: relative tolerance of the GMRES algorithm (default is 1e-10);
    • dealias :: Int: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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.
    source
    WaterWaves1D.SolitaryWhithamType
    SolitaryWhitham(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.

    source
    WaterWaves1D.SolitaryWaveWhithamMethod
    SolitaryWaveWhitham(param; kwargs)

    Compute the Whitham solitary wave with prescribed velocity.

    Argument

    • param :: NamedTuple: parameters of the problem containing velocity c 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 KdV is used);
    • x₀ :: Real: center of solitary wave (if guess is not provided);
    • iterative :: Bool: inverts Jacobian through GMRES if true, LU decomposition if false;
    • verbose :: Bool: prints numerical errors at each step if true;
    • max_iter :: Int: maximum number of iterations of the Newton algorithm;
    • tol :: Real: general tolerance (default is 1e-10);
    • ktol :: Real: tolerance of the Krasny filter (default is 0, i.e. no filtering);
    • gtol :: Real: relative tolerance of the GMRES algorithm (default is 1e-10);
    • dealias :: Int: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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: if true computes the KdV (instead of Whitham) solitary wave.

    Return values

    (u,mesh) with

    • u :: Vector{Float64}: the solution;
    • mesh :: Mesh: mesh collocation points.
    source
    WaterWaves1D.RandomType
    Random(param;args)

    Randomly generated initial data, based on provided (optional arguments) :

    • L is the typical wavelength (default is L=1),
    • s is the (real) Sobolev index regularity (default is s=∞),
    • λ is the length of spatial localization (default is λ=∞, no localization),
    • a is the couple of amplitudes of the surface deformation, and velocity (default is a=(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).

    source
    WaterWaves1D.randomMethod
    random(x;args)

    Randomly generate a vector of size x, based on provided (optional arguments) :

    • L is the typical wavelength (default is L=1),
    • s is the (real) Sobolev index regularity (default is s=∞, smooth data),
    • λ is the length of spatial localization (default is λ=∞, no localization),
    • a is the amplitude of the returned vector (default is a=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.

    source

    Models

      WaterWaves1D.AbstractModelType

      Abstract type whose subtypes are the models from which initial-value problems can be built, through Problem( model :: AbstractModel, initial :: InitialData, param :: NamedTuple )

      source
      WaterWaves1D.AiryType
      Airy(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.
      • 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!:

      1. a function Airy.f! to be called in explicit time-integration solvers;
      2. a function Airy.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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 points x;
        • v is the derivative of the trace of the velocity potential at x.
      source
      WaterWaves1D.AkersNichollsType
      AkersNicholls(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 if IL=true (or μ=Inf, or ν=0), in which case ϵ is the steepness parameter. Default is false.
      • dealias: dealiasing with 1/3 Orlicz rule if true or no dealiasing if false (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!:

      1. a function AkersNicholls.f! to be called in explicit time-integration solvers;
      2. a function AkersNicholls.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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 points x;
        • v is the derivative of the trace of the velocity potential at x.
      source
      WaterWaves1D.BoussinesqType
      Boussinesq(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) and b (default is +1/3) which determine the model solved. You need a+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 is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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!:

      1. a function Boussinesq.f! to be called in explicit time-integration solvers;
      2. a function Boussinesq.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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 points x;
        • v is the derivative of the trace of the velocity potential at x.
      source
      WaterWaves1D.IsobeKakinumaType
      IsobeKakinuma(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 NamedTuples) 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 if true, LU decomposition if false (default is true);
      • precond: use a (left) preconditioner for GMRES if true (default), choose precond as the preconditioner if provided;
      • gtol: relative tolerance of the GMRES algorithm (default is 1e-14);
      • restart: the corresponding option of the GMRES algorithm (default is 100);
      • maxiter: the corresponding option of GMRES (default is nothing);
      • ktol: tolerance of the Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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!:

      1. a function IsobeKakinuma.f! to be called in explicit time-integration solvers;
      2. a function IsobeKakinuma.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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 points x;
        • v is the derivative of the trace of the velocity potential at x.
      source
      WaterWaves1D.MatsunoType
      Matsuno(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 if IL=true (or μ=Inf, or ν=0), in which case ϵ is the steepness parameter. Default is false.
      • mesh: the mesh of collocation points. By default, mesh = Mesh(param);
      • dealias: dealiasing with 1/3 Orlicz rule if true or no dealiasing if false (by default);
      • label: a label for future references (default is "Matsuno");

      Return values

      Generate necessary ingredients for solving an initial-value problem via solve!:

      1. a function Matsuno.f! to be called in explicit time-integration solvers;
      2. a function Matsuno.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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 points x;
        • v is the derivative of the trace of the velocity potential at x.
      source
      WaterWaves1D.NonHydrostaticType
      NonHydrostatic(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 if true, LU decomposition if false (default is true);
      • precond: use a (left) preconditioner for GMRES if true (default), choose precond as the preconditioner if provided;
      • gtol: relative tolerance of the GMRES algorithm (default is 1e-14);
      • restart: the corresponding option of the GMRES algorithm (default is 100);
      • maxiter: the corresponding option of GMRES (default is nothing);
      • ktol: tolerance of the Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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!:

      1. a function NonHydrostatic.f! to be called in explicit time-integration solvers;
      2. a function NonHydrostatic.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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 points x;
        • v is the derivative of the trace of the velocity potential at x;
      4. additionally, a handy function NonHydrostatic.mapfrofull which from data matrix returns the Tuple of real vectors (η,v,u), where
        • u corresponds to the layer-averaged velocity.
      source
      WaterWaves1D.SVNType
      SVN(param;kwargs)

      Define an object of type AbstractModel in view of solving the initial-value problem for a the multilayer Saint-Venant model, possibly with viscosity and/or diffusivity.

      Argument

      param is of type NamedTuple and must contain

      • 'M': the number of layers in the model;
      • dimensionless parameters ϵ (nonlinearity), ν (viscosity) and κ (diffusivity);
      • 'ρ' a Vector of size M, the density of each layer from top to bottom;
      • 'δ' a Vector of size M, the depth of each layer from top to bottom;
      • 'u₀' a Vector of size M, the horizontal velocity of the background shear flow at each layer from top to bottom;
      • 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 is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, i.e. no dealiasing);
      • label: a label for future references (default is "SVN");

      Return values

      Generate necessary ingredients for solving an initial-value problem via solve!:

      1. a function WhithamBoussinesq.f! to be called in explicit time-integration solvers;
      2. a function WhithamBoussinesq.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed.
      3. 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`.
      source
      WaterWaves1D.SaintVenantType
      SaintVenant(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 is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, i.e. no dealiasing);
      • label: a label for future references (default is "Saint-Venant");

      Return values

      Generate necessary ingredients for solving an initial-value problem via solve!:

      1. a function SaintVenant.f! to be called in explicit time-integration solvers;
      2. a function SaintVenant.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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 points x;
        • v is the derivative of the trace of the velocity potential at x.
      source
      WaterWaves1D.SerreGreenNaghdiType
      SerreGreenNaghdi(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 if true, LU decomposition if false (default is true);
      • precond: use a (left) preconditioner for GMRES if true (default), choose precond as the preconditioner if provided;
      • gtol: relative tolerance of the GMRES algorithm (default is 1e-14);
      • restart: the corresponding option of the GMRES algorithm (default is 100);
      • maxiter: the corresponding option of GMRES (default is nothing);
      • ktol: tolerance of the Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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!:

      1. a function SerreGreenNaghdi.f! to be called in explicit time-integration solvers;
      2. a function SerreGreenNaghdi.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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 points x;
        • v is the derivative of the trace of the velocity potential at x;
      4. additionally, a handy function SerreGreenNaghdi.mapfrofull which from data matrix returns the Tuple of real vectors (η,v,u), where
        • u corresponds to the layer-averaged velocity.
      source
      WaterWaves1D.SquareRootDepthType
      SquareRootDepth(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 if true, LU decomposition if false (default is true);
      • precond: use a (left) preconditioner for GMRES if true (default), choose precond as the preconditioner if provided;
      • gtol: relative tolerance of the GMRES algorithm (default is 1e-14);
      • restart: the corresponding option of the GMRES algorithm (default is 100);
      • maxiter: the corresponding option of GMRES (default is nothing);
      • ktol: tolerance of the Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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!:

      1. a function SquareRootDepth.f! to be called in explicit time-integration solvers;
      2. a function SquareRootDepth.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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 points x;
        • v is the derivative of the trace of the velocity potential at x;
      4. additionally, a handy function SquareRootDepth.mapfrofull which from data matrix returns the Tuple of real vectors (η,v,u), where
        • u corresponds to the layer-averaged velocity.
      source
      WaterWaves1D.ToyType
      Toy(param;kwargs)

      Define an object of type AbstractModel which is a toy model for three-scale singular problems

      Argument

      param is of type NamedTuple and must contain

      • the relaxation parameters a and δ;
      • 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);
      • h: the static component function. By default, h = sin;
      • ktol: tolerance of the Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, i.e. no dealiasing);
      • label: a label for future references (default is `"toy");

      Return values

      Generate necessary ingredients for solving an initial-value problem via solve!:

      1. a function toy.f! to be called in explicit time-integration solvers;
      2. a function toy.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. a function toy.mapfro which from such data matrix returns the Tuple of real vectors (h,u,x), where
      - `h` is the values of the static component at `x`;
      - `u` is the values of the real part of solutions at collocation points `x`;
      1. a function toy.mapfrofull which from such data matrix returns the Tuple of real vectors (h,u,x), where
      - `h` is the values of the static component at `x`;
      - `u` is the values of solutions at collocation points `x`;
      source
      WaterWaves1D.WWnType
      WWn(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 if IL=true (or μ=Inf, or ν=0), in which case ϵ is the steepness parameter. Default is false.
      • n :: Int: the order of the expansion; linear system if 1, quadratic if 2, cubic if 3, quartic if 4 (default and other values yield 2);
      • δ and m: parameters of the rectifier operator, set as k->min(1,|δ*k|^m) or k->min(1,|δ*k|^m[1]*exp(1-|δ*k|^m[2])) if m 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 is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, i.e. no dealiasing);
      • label: a label for future references (default is "WWn" with n the order of the expansion);

      Return values

      Generate necessary ingredients for solving an initial-value problem via solve!:

      1. a function WWn.f! to be called in explicit time-integration solvers (also WWn.f1! and WWn.f2! for the symplectic Euler solver);
      2. a function WWn.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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`.
      source
      WaterWaves1D.WaterWavesType
      WaterWaves(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 if IL=true, in which case ϵ is the steepness parameter. Default is false.
      • 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 problem F(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;
      • tol: (relative) tolerance of the fix-point algorithm (default is 1e-16);
      • maxiter: the maximal number of iteration in the fix-point algorithm (default is 100);
      • ktol: tolerance of the low-pass Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, i.e. no dealiasing);
      • label: a label for future references (default is "water waves");
      • verbose: prints information if true (default is true).

      Return values

      Generate necessary ingredients for solving an initial-value problem via solve!:

      1. a function WaterWaves.f! to be called in the explicit time-integration solver (also WaterWaves.f1! and WaterWaves.f2! for the symplectic Euler solver);
      2. a function WaterWaves.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. a function WaterWaves.mapfro which from such data matrix returns the Tuple of real vectors (η,v,x), where
        • x 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`.
      source
      WaterWaves1D.WhithamType
      Whitham(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: if true (default is false), compute the standard KdV equations instead (see KdV(param;kwargs));
      • mesh: the mesh of collocation points. By default, mesh = Mesh(param);
      • ktol: tolerance of the low-pass Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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!:

      1. a function Whitham.f! to be called in explicit time-integration solvers;
      2. a function Whitham.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed.
      3. 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`.
      source
      WaterWaves1D.WhithamBoussinesqType
      WhithamBoussinesq(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: if true (default is false), compute the standard Boussinesq system instead (see Boussinesq(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.
      • mesh: the mesh of collocation points. By default, mesh = Mesh(param);
      • ktol: tolerance of the low-pass Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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!:

      1. a function WhithamBoussinesq.f! to be called in explicit time-integration solvers;
      2. a function WhithamBoussinesq.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed.
      3. 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`.
      source
      WaterWaves1D.WhithamGreenNaghdiType
      WhithamGreenNaghdi(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: if true (default is false), compute the Serre-Green-Naghdi (SGN) instead of Whitham-Green-Naghdi (WGN) system (see SerreGreenNaghdi(param;kwargs));
      • mesh: the mesh of collocation points. By default, mesh = Mesh(param);
      • iterative: solve the elliptic problem through GMRES if true, LU decomposition if false (default is true);
      • precond: use a (left) preconditioner for GMRES if true (default), choose precond as the preconditioner if provided;
      • gtol: relative tolerance of the GMRES algorithm (default is 1e-14);
      • restart: the corresponding option of the GMRES algorithm (default is 100);
      • maxiter: the corresponding option of GMRES (default is nothing);
      • ktol: tolerance of the Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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!:

      1. a function WhithamGreenNaghdi.f! to be called in explicit time-integration solvers;
      2. a function WhithamGreenNaghdi.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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`;
      1. additionally, a handy function WhithamGreenNaghdi.mapfrofull which from data matrix returns the Tuple of real vectors (η,v,u), where
        • u corresponds to the layer-averaged velocity.
      source
      WaterWaves1D.modifiedMatsunoType
      modifiedMatsuno(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 if IL=true (or μ=Inf, or ν=0), in which case ϵ is the steepness parameter. Default is false.
      • mesh: the mesh of collocation points. By default, mesh = Mesh(param);
      • ktol: tolerance of the low-pass Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, 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!:

      1. a function modifiedMatsuno.f! to be called in explicit time-integration solvers;
      2. a function modifiedMatsuno.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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 points x;
        • v is the derivative of the trace of the velocity potential at x.
      source
      WaterWaves1D.relaxedGreenNaghdiType
      relaxedGreenNaghdi(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: if true (default is false), 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 is 1);
      • mesh: the mesh of collocation points. By default, mesh = Mesh(param);
      • iterative: solve the elliptic problem (to construct initial data) through GMRES if true, LU decomposition if false (default is true);
      • precond: use a (left) preconditioner for GMRES if true (default), choose precond as the preconditioner if provided;
      • gtol: relative tolerance of the GMRES algorithm (default is 1e-14);
      • restart: the corresponding option of the GMRES algorithm (default is 100);
      • maxiter: the corresponding option of GMRES (default is nothing);
      • ktol: tolerance of the Krasny filter (default is 0, i.e. no filtering);
      • dealias: dealiasing with Orlicz rule 1-dealias/(dealias+2) (default is 0, i.e. no dealiasing);
      • label: a label for future references (default is "Favrie-Gavrilyuk" if FG==true, "Escalante-Dumbser-Castro" otherwise);

      Return values

      Generate necessary ingredients for solving an initial-value problem via solve!:

      1. a function relaxedGreenNaghdi.f! to be called in explicit time-integration solvers;
      2. a function relaxedGreenNaghdi.mapto which from (η,v) of type InitialData provides the raw data matrix on which computations are to be executed;
      3. 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`;
      1. 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.
      source

      Solvers

        WaterWaves1D.EulerType
        Euler(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

        1. an object of type AbstractModel;
        2. an Array of size (N,datasize) where N is the number of collocation points and datasize the number of equations solved;
        3. (param,datasize) where param is aNamedTuplecontaining a keyN, anddatasizea 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..

        source
        WaterWaves1D.EulerSympType
        EulerSymp(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

        1. an object of type AbstractModel;
        2. an Array of size (N,2) where N is the number of collocation points;
        3. a NamedTuple containing a key N.

        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..

        source
        WaterWaves1D.Euler_naiveType
        Euler_naive()

        Runge-Kutta fourth order solver.

        A naive version of Euler, without argument since no pre-allocation is performed.

        source
        WaterWaves1D.RK4Type
        RK4(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

        1. an object of type AbstractModel, typically the model you will solve with the solver;
        2. an Array of size (N,datasize) where N is the number of collocation points and datasize the number of equations solved;
        3. (param,datasize) where param is a NamedTuple containing a key N, and datasize a 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..

        source
        WaterWaves1D.RK4_naiveType
        RK4_naive()

        Runge-Kutta fourth order solver.

        A naive version of RK4, without argument since no pre-allocation is performed.

        source

        Structures

        WaterWaves1D.InitType
        Init(data ; fast, label)

        Generate an initial data to be used in the function Problem.

        data should contain either

        • a function η and a function v (in this order);
        • a mesh (of type Mesh) and two vectors representing η(mesh.x) and v(mesh.x) (in this order);
        • an array of collocation points and two vectors representing η(x) and v(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".

        source
        WaterWaves1D.InitialDataType

        Abstract type defining initial data from which initial-value problems can be built, through Problem( model :: AbstractModel, initial :: InitialData, param :: NamedTuple )

        source
        WaterWaves1D.MeshType
        Mesh(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 containing N and xmin and xmax; or
        • a NamedTuple containing N and L (in which case xmin=-L and xmax=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.
        source
        WaterWaves1D.ProblemType
        Problem( 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 function Times.

        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 is true).
        source
        WaterWaves1D.solve!Method
        solve!( 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).

        source
        WaterWaves1D.solve!Method
        solve!( 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).

        source
        WaterWaves1D.DataType
        Data( 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 matrix m;
        • (data.datalength,data.datasize)=size(m) where datalength is the number of computed modes, and datasize the number of involved equations, typically 2.
        source
        WaterWaves1D.TimesType
        Times(param; ns, Ns)

        Constructs a mesh of times, to be used in initial-value problems (see Problem).

        Arguments

        param is either

        • a NamedTuple containing dt the timestep and T the final time of comuptation; or
        • a vector of computed times.

        Optional keyword arguments

        • ns : data are stored every ns 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.
        source

        Tools

        WaterWaves1D.energyMethod
        energy(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.

        source
        WaterWaves1D.energydiffMethod
        energydiff(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).

        source
        WaterWaves1D.interpolateMethod
        interpolate(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.

        source
        WaterWaves1D.interpolateMethod
        interpolate(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.

        source
        WaterWaves1D.massMethod
        mass(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.

        source
        WaterWaves1D.massdiffMethod
        massdiff(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).

        source
        WaterWaves1D.momentumMethod
        momentum(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.

        source
        WaterWaves1D.momentumdiffMethod
        momentumdiff(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).

        source
        WaterWaves1D.solutionMethod
        solution(pb::Problem;T,x,interpolation)

        Give the solution of a solved initial-value problem at a given time T.

        Arguments

        • Argument pb is of type Problem.
        • 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 vector x.
        • Keyword argument interpolation is optional, if an integer is provided the solution is interpolated on as many collocation points (if true, then the default value 2^3 is chosen).
        • Keyword argument raw is optional, if set to true then (U,t) with U the raw data and t the time is returned (default is false).

        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 provided T).
        source

        Load and save

        Base.dumpMethod
        dump(file_name :: String, data :: Data)

        Save data to the file with name file_name (and extension ".h5").

        source
        Base.dumpMethod
        dump(file_name  :: String, problem :: Problem)

        Save the data of problem to the file with name file_name (and extension ".h5").

        source
        Base.dumpMethod
        dump(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").

        source
        WaterWaves1D.load_data!Method
        load_data!(file_name :: String, problem :: Problem)

        Fills problem with raw data extracted from the file with name file_name (and extension ".h5").

        source
        WaterWaves1D.load_initMethod
        load_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.

        source