2D Heat Equation
Overview
This is the 2D version of 1D heat equation. It is defined by
\[u_t = \mu (u_{xx} + u_{yy})\]
where $x,y\in[0,L]$ is the temperature and $\mu$ is the thermal diffusivity parameter.
Finite Difference Model
The finite difference approach is similar to the 1D version but with additional complications due to the addition of an extra dimension. If the 2D domain is spatially discretized using $N$ and $M$ grid points for the $x$ and $y$ directions, respectively, then the toeplitz matrices corresponding to each $x$ and $y$ directions are identical to that of the 1D case which is defined by
\[\mathbf{A}_x\in\mathbb{R}^{N\times N} \qquad \text{and} \qquad \mathbf{A}_y\in\mathbb{R}^{M\times M}.\]
However, to construct the matrix for the overall system we utilize the Kronecker product. Define the state vector $\mathbf{z}$ which flattens the 2D grid into a vector then the linear matrix becomes
\[\mathbf{A} = \mathbf{A}_y \otimes \mathbf{I}_{N} + \mathbf{I}_M \otimes \mathbf{A}_x\]
and the $\mathbf{B}$ matrix will be constructed such that they add the inputs to the appropriate indices of the flattened state vector $\mathbf z$.
Thus, we arrive at a $N$-dimensional linear time-invariant (LTI) system:
\[\dot{\mathbf{u}}(t) = \mathbf{A}\mathbf{u}(t) + \mathbf{B}\mathbf{w}(t)\]
We then consider the numerical integration scheme. For our numerical integration we can consider three approaches
- Forward Euler
- Backward Euler
- Crank-Nicolson
Refer to 1D heat equation for details on the numerical integration schemes.
Example
using CairoMakie
using LinearAlgebra
using PolynomialModelReductionDataset: Heat2DModel
using UniqueKronecker: invec
# Setup
Ω = ((0.0, 1.0), (0.0, 1.0))
Nx = 2^5
Ny = 2^5
heat2d = Heat2DModel(
spatial_domain=Ω, time_domain=(0,2),
Δx=(Ω[1][2] + 1/Nx)/Nx, Δy=(Ω[2][2] + 1/Ny)/Ny, Δt=1e-3,
diffusion_coeffs=0.1, BC=(:dirichlet, :dirichlet)
)
xgrid0 = heat2d.yspan' .* ones(heat2d.spatial_dim[1])
ygrid0 = ones(heat2d.spatial_dim[2])' .* heat2d.xspan
ux0 = sin.(2π * xgrid0) .* cos.(2π * ygrid0)
heat2d.IC = vec(ux0) # initial condition
# Boundary condition
Ubc = [1.0, 1.0, -1.0, -1.0]
Ubc = repeat(Ubc, 1, heat2d.time_dim)
# Operators
A, B = heat2d.finite_diff_model(heat2d, heat2d.diffusion_coeffs)
# Integrate
U = heat2d.integrate_model(
heat2d.tspan, heat2d.IC, Ubc;
operators=[A,B], system_input=true, integrator_type=:BackwardEuler
)
U2d = invec.(eachcol(U), heat2d.spatial_dim...)
fig1 = Figure()
ax1 = Axis3(fig1[1, 1], xlabel="x", ylabel="y", zlabel="u(x,y,t)")
sf = surface!(ax1, heat2d.xspan, heat2d.yspan, U2d[1])
Colorbar(fig1[1, 2], sf)
fig1

fig2 = Figure()
ax2 = Axis3(fig2[1, 1], xlabel="x", ylabel="y", zlabel="u(x,y,t)")
sf = surface!(ax2, heat2d.xspan, heat2d.yspan, U2d[end])
Colorbar(fig2[1, 2], sf)
fig2

API
PolynomialModelReductionDataset.Heat2D.Heat2DModel
— Typemutable struct Heat2DModel <: AbstractModel
2 Dimensional Heat Equation Model
\[\frac{\partial u}{\partial t} = \mu\left(\frac{\partial^2 u}{\partial x^2} + \frac{\partial^2 u}{\partial y^2}\right)\]
Fields
spatial_domain::Tuple{Tuple{<:Real,<:Real}, Tuple{<:Real,<:Real}}
: spatial domain (x, y)time_domain::Tuple{Real,Real}
: temporal domainparam_domain::Tuple{Real,Real}
: parameter domainΔx::Real
: spatial grid size (x-axis)Δy::Real
: spatial grid size (y-axis)Δt::Real
: temporal step sizespatial_dim::Tuple{Int64,Int64}
: spatial dimension x and ytime_dim::Int64
: temporal dimensionparam_dim::Int64
: parameter dimensionBC::Tuple{Symbol,Symbol}
: boundary conditionIC::Array{Float64}
: initial conditiondiffusion_coeffs::Union{AbstractArray{<:Real},Real}
: diffusion coefficientsxspan::Vector{Float64}
: spatial grid points (x-axis)yspan::Vector{Float64}
: spatial grid points (y-axis)tspan::Vector{Float64}
: temporal pointsfinite_diff_model::Function
: finite difference modelintegrate_model::Function
: integrate model
PolynomialModelReductionDataset.Heat2D
— Module2D Heat Equation Model
PolynomialModelReductionDataset.Heat2D.finite_diff_model
— Methodfinite_diff_model(model, μ)
Generate A and B matrices for the 2D heat equation.
Arguments
model::heat2d
: 2D heat equation modelμ::Real
: diffusion coefficient
Returns
A::Matrix{Float64}
: A matrixB::Matrix{Float64}
: B matrix
PolynomialModelReductionDataset.Heat2D.integrate_model
— Methodintegrate_model(A, B, U, tdata, IC)
Integrate the 2D heat equation model.
Arguments
A::Matrix{Float64}
: A matrixB::Matrix{Float64}
: B matrixU::Vector{Float64}
: input vectortdata::Vector{Float64}
: time pointsIC::Vector{Float64}
: initial condition
Returns
state::Matrix{Float64}
: state matrix
PolynomialModelReductionDataset.Heat2D.integrate_model
— Methodintegrate_model(tdata, u0)
integrate_model(tdata, u0, input; kwargs...)
Integrate the 1D Heat Equation Model using 3 different methods:
- Forward Euler
- Backward Euler
- Crank Nicolson
Arguments
tdata::AbstractVector{T}
: time datax0::AbstractVector{T}
: initial conditionu::AbstractArray{T}=[]
: input data
Keyword Arguments
operators
: operators A and Bsystem_input::Bool=false
: system input flagintegrator_type::Symbol=:ForwardEuler
: integrator type
Returns
x::Array{T,2}
: integrated model states
Notes
- Input is assumed to be a matrix of size (spatial dimension x time dimension). You will receive a warning if the input is a tall vector/column vector.
operators
should be in the order of [A, B] ifsystem_input
is true.