pymgrit.heat package¶
Submodules¶
pymgrit.heat.heat_1d module¶
Vector and application class for the 1D heat equation
- class pymgrit.heat.heat_1d.Heat1D(*args, **kwargs)¶
Bases:
pymgrit.core.application.Application
- Application class for the heat equation in 1D space,
u_t - a*u_xx = b(x,t), a > 0, x in [x_start,x_end], t in [0,T],
with homogeneous Dirichlet boundary conditions in space.
- compute_matrix()¶
Define spatial discretization matrix for 1D heat equation
- Second-order central finite differences with matrix stencil
(a / dx^2) * [-1 2 -1]
- step(u_start: pymgrit.heat.heat_1d.VectorHeat1D, t_start: float, t_stop: float) pymgrit.heat.heat_1d.VectorHeat1D ¶
- Time integration routine for 1D heat equation example problem:
Backward Euler (BDF1)
- One-step method
u_i = (I + dt*L)^{-1} * (u_{i-1} + dt*b_i),
where L = self.space_disc is the spatial discretization operator
- Parameters
u_start – approximate solution for the input time t_start
t_start – time associated with the input approximate solution u_start
t_stop – time to evolve the input approximate solution to
- Returns
approximate solution at input time t_stop
- class pymgrit.heat.heat_1d.VectorHeat1D(size)¶
Bases:
pymgrit.core.vector.Vector
Vector class for the 1D heat equation
- clone()¶
Clone vector object
- Return type
vector object with zero values
- clone_rand()¶
Initialize vector object with random values
- Return type
vector object with random values
- clone_zero()¶
Initialize vector object with zeros
- Return type
vector object with zero values
- get_values()¶
Get vector data
- Returns
values of vector object
- norm()¶
Norm of a vector object
- Returns
2-norm of vector object
- pack()¶
Pack data
- Returns
values of vector object
- set_values(values)¶
Set vector data
- Parameters
values – values for vector object
- unpack(values)¶
Unpack and set data
- Parameters
values – values for vector object
pymgrit.heat.heat_1d_2pts_bdf1 module¶
Application class for 1D heat problem using BDF1 time integration
Note: values at two consecutive time points are grouped as pairs
- class pymgrit.heat.heat_1d_2pts_bdf1.Heat1DBDF1(*args, **kwargs)¶
Bases:
pymgrit.core.application.Application
- Application class for the heat equation in 1D space,
u_t - a*u_xx = b(x,t), a > 0, x in [x_start,x_end], t in [0,T],
with homogeneous Dirichlet boundary conditions in space
- compute_matrix()¶
Define spatial discretization matrix for heat equation problem.
- Discretization is centered finite differences with matrix stencil
(a / dx^2) * [-1 2 -1]
- step(u_start: pymgrit.heat.vector_heat_1d_2pts.VectorHeat1D2Pts, t_start: float, t_stop: float) pymgrit.heat.vector_heat_1d_2pts.VectorHeat1D2Pts ¶
- Time integration routine for 1D heat equation example problem:
Backward Euler (BDF1)
- One-step method
u_i = (I + dt*L)^{-1} * (u_{i-1} + dt*b_i),
where L = self.space_disc is the spatial discretization operator
- Note: step takes two BDF1 steps
one step from (t_start + dtau) to t_stop
one step from t_stop to (t_stop + dtau)
- Parameters
u_start – approximate solution for the input time t_start
t_start – time associated with the input approximate solution u_start
t_stop – time to evolve the input approximate solution to
- Returns
approximate solution at input time t_stop
pymgrit.heat.heat_1d_2pts_bdf2 module¶
Application class for 1D heat problem using BDF2 time integration
Note: values at two consecutive time points are grouped as pairs
- class pymgrit.heat.heat_1d_2pts_bdf2.Heat1DBDF2(*args, **kwargs)¶
Bases:
pymgrit.core.application.Application
- Application class for the heat equation in 1D space,
u_t - a*u_xx = b(x,t), a > 0, x in [x_start, x_end], t in [t_start, t_end],
with homogeneous Dirichlet boundary conditions in space
- compute_matrix()¶
Define spatial discretization matrix for heat equation problem.
- Discretization is centered finite differences with matrix stencil
(a / dx^2) * [-1 2 -1]
- step(u_start: pymgrit.heat.vector_heat_1d_2pts.VectorHeat1D2Pts, t_start: float, t_stop: float) pymgrit.heat.vector_heat_1d_2pts.VectorHeat1D2Pts ¶
- Time integration routine for 1D heat equation:
BDF2
Two-step method on variably spaced grid with spacing tau_i = t_i - t_{i-1}. In time-based stencil notation, we have at time point t_i
[r_i^2/(tau_i*(1+r_i))*I, -((1+r_i)/tau_i)*I, (1+2r_i)/(tau_i*(1+r_i))*I + L, 0, 0],
where L = self.space_disc is the spatial discretization operator and r_i = tau_i/tau_{i-1}
- Note: For the pair associated with input time t_stop
update at t_stop involves values at t_start and (t_start + dtau)
update at t_stop + dtau involves values at (t_start + dtau) and t_stop
- Parameters
u_start – approximate solution for the input time t_start
t_start – time associated with the input approximate solution u_start
t_stop – time to evolve the input approximate solution to
- Returns
approximate solution for the input time t_stop
pymgrit.heat.heat_2d module¶
Vector and application class for the 2D heat equation
- Time integration methods:
Forward Euler Backward Euler Crank-Nicolson
- class pymgrit.heat.heat_2d.Heat2D(*args, **kwargs)¶
Bases:
pymgrit.core.application.Application
- Application class for the heat equation in 2D space,
u_t - a(u_xx + u_yy) = b(x,y,t), a > 0, in [x_start, x_end] x [y_start, y_end] x (t_start, t_end],
with homogeneous Dirichlet boundary conditions in space.
- compute_matrix()¶
Define spatial discretization matrix for 2D heat equation
- Second-order central finite differences with matrix stencil
[ -f_y ] [-f_x 2(f_x + f_y) -f_x] [ -f_y ]
with f_x = (a / dx^2) and f_y = (a / dy^2)
- compute_rhs(u_start: numpy.ndarray, t_start: float, t_stop: float)¶
Right-hand side of spatial system for implicit time integration methods
- Parameters
u_start – approximate solution for the input time t_start
t_start – time associated with the input approximate solution u_start
t_stop – time to evolve the input approximate solution to
- Returns
right-hand side of spatial system at each time step in case of implicit time integration
- step(u_start: pymgrit.heat.heat_2d.VectorHeat2D, t_start: float, t_stop: float) pymgrit.heat.heat_2d.VectorHeat2D ¶
- Time integration routine for 2D heat equation example problem:
Backward Euler (BE) Forward Euler (FE) Crank-Nicolson (CN)
- One-step method
(u_i - u_{i-1})/dt + (theta*L*u_i + (1-theta)*L*u_{i-1} = theta*b_i + (1-theta)*b_{i-1},
- where L = self.space_disc is the spatial discretization operator and
theta = 0 -> FE theta = 1/2 -> CN theta = 1 -> BE
- Parameters
u_start – approximate solution for the input time t_start
t_start – time associated with the input approximate solution u_start
t_stop – time to evolve the input approximate solution to
- Returns
approximate solution at input time t_stop
- class pymgrit.heat.heat_2d.VectorHeat2D(nx, ny)¶
Bases:
pymgrit.core.vector.Vector
Vector class for the 2D heat equation
- clone()¶
Clone vector object
- Return type
vector object with zero values
- clone_rand()¶
Initialize vector object with random values
- Return type
vector object with random values
- clone_zero()¶
Initialize vector object with zeros
- Return type
vector object with zero values
- get_values()¶
Get vector data
- Returns
values of vector object
- norm()¶
Norm of a vector object
- Returns
2-norm of vector object
- pack()¶
Pack data
- Returns
values of vector object
- set_values(values)¶
Set vector data
- Parameters
values – values for vector object
- unpack(values)¶
Unpack and set data
- Parameters
values – values for vector object
pymgrit.heat.vector_heat_1d_2pts module¶
Vector class for 1D heat problem
Note: values at two consecutive time points are grouped as pairs
- class pymgrit.heat.vector_heat_1d_2pts.VectorHeat1D2Pts(size, dtau)¶
Bases:
pymgrit.core.vector.Vector
Vector class for grouping values at two consecutive time points
- clone()¶
Clone vector object
- Returns
vector object with zero values
- clone_rand()¶
Initialize vector object with random values
- Returns
vector object with random values
- clone_zero()¶
Initialize vector object with zeros
- Returns
vector object with zero values
- get_values()¶
Get vector data
- Returns
tuple of values of member variables
- norm()¶
Norm of a vector object
- Returns
2-norm of vector object
- pack()¶
Pack data
- Returns
values of vector object
- set_values(first_time_point, second_time_point, dtau)¶
Set vector data
- Parameters
first_time_point – values for first time point
second_time_point – values for second time point
dtau – time-step size within pair
- unpack(values)¶
Unpack and set data
- Parameters
values – values for vector object