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
Example block output
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
Example block output

API

PolynomialModelReductionDataset.Heat2D.Heat2DModelType
mutable 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 domain
  • param_domain::Tuple{Real,Real}: parameter domain
  • Δx::Real: spatial grid size (x-axis)
  • Δy::Real: spatial grid size (y-axis)
  • Δt::Real: temporal step size
  • spatial_dim::Tuple{Int64,Int64}: spatial dimension x and y
  • time_dim::Int64: temporal dimension
  • param_dim::Int64: parameter dimension
  • BC::Tuple{Symbol,Symbol}: boundary condition
  • IC::Array{Float64}: initial condition
  • diffusion_coeffs::Union{AbstractArray{<:Real},Real}: diffusion coefficients
  • xspan::Vector{Float64}: spatial grid points (x-axis)
  • yspan::Vector{Float64}: spatial grid points (y-axis)
  • tspan::Vector{Float64}: temporal points
  • finite_diff_model::Function: finite difference model
  • integrate_model::Function: integrate model
source
PolynomialModelReductionDataset.Heat2D.integrate_modelMethod
integrate_model(A, B, U, tdata, IC)

Integrate the 2D heat equation model.

Arguments

  • A::Matrix{Float64}: A matrix
  • B::Matrix{Float64}: B matrix
  • U::Vector{Float64}: input vector
  • tdata::Vector{Float64}: time points
  • IC::Vector{Float64}: initial condition

Returns

  • state::Matrix{Float64}: state matrix
source
PolynomialModelReductionDataset.Heat2D.integrate_modelMethod
integrate_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 data
  • x0::AbstractVector{T}: initial condition
  • u::AbstractArray{T}=[]: input data

Keyword Arguments

  • operators: operators A and B
  • system_input::Bool=false: system input flag
  • integrator_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] if system_input is true.
source