Damped Gardner-Burgers Equation
Overview
This equation is a modified version of the Gardner equation where the diffusion term $u_{xx}$ and damping term $u$ is added for stability. The equation is defined by
\[u_t = -\alpha u_{xxx} + \beta uu_x + \gamma u^2 u_{x} + \delta u_{xx} + \epsilon u\]
Model
The damped Gardner-Burgers equation becomes the exact same as the original Gardner equation:
\[\dot{\mathbf{u}}(t) = \mathbf{Au}(t) + \mathbf{H}(\mathbf{u}(t) \otimes \mathbf{u}(t)) + \mathbf{G}(\mathbf{u}(t) \otimes \mathbf{u}(t) \otimes \mathbf{u}(t)) + \mathbf{Bw}(t)\]
or
\[\dot{\mathbf{u}}(t) = \mathbf{Au}(t) + \mathbf{F}(\mathbf{u}(t) \oslash \mathbf{u}(t)) + \mathbf{E}(\mathbf{u}(t) \oslash \mathbf{u}(t) \oslash \mathbf{u}(t)) + \mathbf{Bw}(t)\]
where
- $\mathbf{u}\in\mathbb{R}^N$: the state vector
- $\mathbf{w}\in\mathbb{R}^m$: the input vector
- $\mathbf{A}\in\mathbb{R}^{N\times N}$: the linear state matrix
- $\mathbf{H}\in\mathbb{R}^{N\times N^2}$: the quadratic state matrix with redundancy
- $\mathbf{F}\in\mathbb{R}^{N\times N(N+1)/2}$: the quadratic state matrix without redundancy
- $\mathbf{G}\in\mathbb{R}^{N\times N^3}$: the cubic state matrix with redundancy
- $\mathbf{E}\in\mathbb{R}^{N\times N(N+1)(N+2)/6}$: the cubic state matrix without redundancy
- $\mathbf{B}\in\mathbb{R}^{N\times m}$: the control input matrix
Numerical Integration
For the numerical integration we consider two methods:
- Semi-Implicit Euler (SIE)
- Crank-Nicolson Adam-Bashforth (CNAB)
For the exact expressions of the time-stepping check Allen-Cahn equation.
Examples
using CairoMakie
using LinearAlgebra
using PolynomialModelReductionDataset: DampedGardnerBurgersModel
# Setup
Ω = (0.0, 3.0)
Nx = 2^8; dt = 1e-3
dgb = DampedGardnerBurgersModel(
spatial_domain=Ω, time_domain=(0.0, 3.0), Δx=(Ω[2] + 1/Nx)/Nx, Δt=dt,
params=Dict(:a => 1, :b => 3, :c => 5, :d => 0.2, :e => 0.5), BC=:dirichlet,
)
DS = 100
dgb.IC = 2 * cos.(2π * dgb.xspan / (Ω[2] - Ω[1])) # + 0.5 * cos.(4π * dgb.xspan / (Ω[2] - Ω[1]))
Ubc1 = 0.5ones(1,dgb.time_dim)
Ubc2 = -0.5ones(1,dgb.time_dim)
Ubc = [Ubc1; Ubc2]
# Operators
A, F, E, B = dgb.finite_diff_model(dgb, dgb.params)
# Integrate
U = dgb.integrate_model(
dgb.tspan, dgb.IC, Ubc;
linear_matrix=A, quadratic_matrix=F, cubic_matrix=E, control_matrix=B,
system_input=true, integrator_type=:CNAB,
)
# Surface plot
fig3, _, sf = CairoMakie.surface(dgb.xspan, dgb.tspan[1:DS:end], U[:, 1:DS:end],
axis=(type=Axis3, xlabel=L"x", ylabel=L"t", zlabel=L"u(x,t)"))
CairoMakie.Colorbar(fig3[1, 2], sf)
fig3

# Flow field
fig4, ax, hm = CairoMakie.heatmap(dgb.xspan, dgb.tspan[1:DS:end], U[:, 1:DS:end])
ax.xlabel = L"x"
ax.ylabel = L"t"
CairoMakie.Colorbar(fig4[1, 2], hm)
fig4

API
PolynomialModelReductionDataset.DampedGardnerBurgers.DampedGardnerBurgersModel
— Typemutable struct DampedGardnerBurgersModel <: AbstractModel
Damped Gardner-Burgers equation model
\[\frac{\partial u}{\partial t} = -\alpha\frac{\partial^3 u}{\partial x^3} - \beta u\frac{\partial u}{\partial x} - \gamma u^2\frac{\partial u}{\partial x} + \delta \frac{\partial^2 u}{\partial x^2} + \epsilon u\]
Fields
spatial_domain::Tuple{Real,Real}
: spatial domaintime_domain::Tuple{Real,Real}
: temporal domainparam_domain::Tuple{Real,Real}
: parameter domainΔx::Real
: spatial grid sizeΔt::Real
: temporal step sizeBC::Symbol
: boundary conditionIC::Array{Float64}
: initial conditionxspan::Vector{Float64}
: spatial grid pointstspan::Vector{Float64}
: temporal pointsspatial_dim::Int64
: spatial dimensiontime_dim::Int64
: temporal dimensionparams::Union{Real,AbstractArray{<:Real}}
: parameter vectorparam_dim::Int64
: parameter dimensionfinite_diff_model::Function
: model using Finite Differenceintegrate_model::Function
: model integration
PolynomialModelReductionDataset.DampedGardnerBurgers
— ModuleDamped Gardner-Burgers (DGB) PDE model
PolynomialModelReductionDataset.DampedGardnerBurgers.finite_diff_dirichlet_model
— Methodfinite_diff_dirichlet_model(N, Δx, Δt, params)
Generate A, F, E, B matrices for the Gardner equation for Dirichlet boundary condition.
PolynomialModelReductionDataset.DampedGardnerBurgers.finite_diff_model
— Methodfinite_diff_model(model, params)
Finite Difference Model for Gardner equation
Arguments
model::DampedGardnerBurgersModel
: Gardner modelparams::Real
: params including a, b, c
Returns
- operators
PolynomialModelReductionDataset.DampedGardnerBurgers.finite_diff_periodic_model
— Methodfinite_diff_periodic_model(N, Δx, params)
Generate A, F, E matrices for the Gardner equation for periodic boundary condition (Non-conservative).
PolynomialModelReductionDataset.DampedGardnerBurgers.integrate_model
— Methodintegrate_model(tdata, u0)
integrate_model(tdata, u0, input; kwargs...)
Integrate the Damped Gardner-Burgers model using 2 different methods:
- Semi-Implicit Euler (SIE)
- Crank-Nicolson Adam-Bashforth (CNAB)
Arguments
tdata::AbstractArray{T}
: time datau0::AbstractArray{T}
: initial conditioninput::AbstractArray{T}=[]
: input data
Keyword Arguments
linear_matrix::AbstractArray{T,2}
: linear matrixquadratic_matrix::AbstractArray{T,2}
: quadratic matrixcubic_matrix::AbstractArray{T,2}
: cubic matrixcontrol_matrix::AbstractArray{T,2}
: control matrixsystem_input::Bool=false
: system inputconst_stepsize::Bool=false
: constant step sizeu2_jm1::Union{Nothing,AbstractArray{T}}=nothing
: u2 at j-1u3_jm1::Union{Nothing,AbstractArray{T}}=nothing
: u3 at j-1integrator_type::Symbol=:CNAB
: integrator type
Returns
u::AbstractArray{T}
: solution
Notes
- If
system_input=true
, then the control matrix must be provided - If
const_stepsize=true
, then the time step size is assumed to be constant - If
u2_jm1
andu3_jm1
are provided, then the CNAB scheme is used - If
integrator_type=:SIE
, then the Semi-Implicit Euler scheme is used - If
integrator_type=:CNAB
, then the Crank-Nicolson Adam-Bashforth scheme is used
PolynomialModelReductionDataset.DampedGardnerBurgers.integrate_model_with_control_CNAB
— Methodintegrate_model_with_control_CNAB(
tdata,
u0,
input;
linear_matrix,
quadratic_matrix,
cubic_matrix,
control_matrix,
const_stepsize,
u2_jm1,
u3_jm1
)
Crank-Nicolson Adam-Bashforth (CNAB) scheme with control
PolynomialModelReductionDataset.DampedGardnerBurgers.integrate_model_with_control_SIE
— Methodintegrate_model_with_control_SIE(
tdata,
u0,
input;
linear_matrix,
quadratic_matrix,
cubic_matrix,
control_matrix
)
Semi-Implicit Euler scheme with control
PolynomialModelReductionDataset.DampedGardnerBurgers.integrate_model_without_control_CNAB
— Methodintegrate_model_without_control_CNAB(
tdata,
u0;
linear_matrix,
quadratic_matrix,
cubic_matrix,
const_stepsize,
u2_jm1,
u3_jm1
)
Crank-Nicolson Adam-Bashforth (CNAB) scheme without control
PolynomialModelReductionDataset.DampedGardnerBurgers.integrate_model_without_control_SIE
— Methodintegrate_model_without_control_SIE(
tdata,
u0;
linear_matrix,
quadratic_matrix,
cubic_matrix
)
Semi-Implicit Euler scheme without control