pymgrit package

Subpackages

Module contents

class pymgrit.Advection1D(*args, **kwargs)

Bases: pymgrit.core.application.Application

Application class for the advection problem in 1D space,

u_t + c * u_t = 0,

subject to periodic boundary conditions in space and initial condition u(x, 0) = exp(-x^2).

compute_matrix()

Define spatial discretization matrix for advection problem.

Discretization is first-order upwind in space.

initialise()
Set the initial condition of the 1D advection problem

u(x,0) = exp(-x^2)

step(u_start: pymgrit.advection.advection_1d.VectorAdvection1D, t_start: float, t_stop: float)
Time integration routine for 1D advection problem:

Backward Euler

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

class pymgrit.Application(*args, **kwargs)

Bases: object

Abstract application class for user-defined application classes. Every user-defined application class must inherit from this class.

Required attributes:
  • vector_template

  • vector_t_start

Required functions:
  • step

required_attributes = ['vector_template', 'vector_t_start']
abstract step(u_start: object, t_start: float, t_stop: float) object

Time integration routine for the application

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

property vector_t_start: pymgrit.core.vector.Vector

Getter for attribute vector_t_start

property vector_template: pymgrit.core.vector.Vector

Getter for attribute vector_template

class pymgrit.ArenstorfOrbit(*args, **kwargs)

Bases: pymgrit.core.application.Application

Application for Arenstorf orbit problem,

x’’ = x + 2y’ - b*(x + a)/D_1 - a*(x - b)/D_2, y’’ = y - 2x’ - b*y/D_1 - a*y/D_2

with a = 0.012277471, b = 1 - a,

D_1 = ((x + a)^2 + y^2)^(3/2), D_2 = ((x - a)^2 + y^2)^(3/2)

and ICs

x(0) = 0.994, x’(0) = 0, y(0) = 0, y’(0) = -2.00158510637908

step(u_start: pymgrit.arenstorf_orbit.arenstorf_orbit.VectorArenstorfOrbit, t_start: float, t_stop: float) pymgrit.arenstorf_orbit.arenstorf_orbit.VectorArenstorfOrbit

Time integration routine for the application

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.Brusselator(*args, **kwargs)

Bases: pymgrit.core.application.Application

Application class for Brusselator system,

x’ = A + x^2y - (B + 1)x, y’ = Bx - x^2y,

with A = 1, B = 3, and ICs

x(0) = 0, y(0) = 1

step(u_start: pymgrit.brusselator.brusselator.VectorBrusselator, t_start: float, t_stop: float) pymgrit.brusselator.brusselator.VectorBrusselator

Time integration routine for Brusselator system: RK4

0 |

1 / 2 | 1 / 2 1 / 2 | 0 1 / 2

1 | 0 0 1

——+—————————-
1 / 6 1 / 3 1 / 3 1 / 6
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

class pymgrit.Dahlquist(*args, **kwargs)

Bases: pymgrit.core.application.Application

Application class for Dahlquist’s test problem,

u’ = lambda u,

with lambda = -1 and IC u(0) = 1

step(u_start: pymgrit.dahlquist.dahlquist.VectorDahlquist, t_start: float, t_stop: float) pymgrit.dahlquist.dahlquist.VectorDahlquist
Time integration routine for Dahlquist’s test problem:

BE: Backward Euler FE: Forward Euler TR: Trapezoidal rule MR: implicit Mid-point Rule

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

class pymgrit.GridTransfer

Bases: abc.ABC

Abstract grid transfer class for user-defined grid transfer classes. Every user-defined grid transfer class must inherit from this class. Allows for additional spatial coarsening and refinement between time levels.

Required functions:
  • restriction

  • interpolation

abstract interpolation(u: pymgrit.core.vector.Vector) pymgrit.core.vector.Vector

Spatial interpolation. Receives the spatial solution at a point in time of one level and interpolates it to the same point in time on the next finer level.

Parameters

u – Approximate solution at a time point of one level

Returns

Interpolated input vector u on next finer time grid

Return type

Vector

abstract restriction(u: pymgrit.core.vector.Vector) pymgrit.core.vector.Vector

Spatial restriction. Receives the spatial solution at a point in time of one level and restricts it to the same point in time on the next coarser level.

Parameters

u – Approximate solution at a time point of one level

Returns

Restricted input vector u on next coarser time grid

Return type

Vector

class pymgrit.GridTransferCopy

Bases: pymgrit.core.grid_transfer.GridTransfer

Standard grid transfer. Copies the spatial solution from one level to another. This function is called for every time point.

interpolation(u: pymgrit.core.vector.Vector) pymgrit.core.vector.Vector

Copies the spatial solution at one point in time of one level and restricts it to the same point in time on the next finer level.

Parameters

u – Approximate solution at a time point of one level

Returns

Input vector u at the same time point on next finer time grid

Return type

Vector

restriction(u: pymgrit.core.vector.Vector) pymgrit.core.vector.Vector

Copies the spatial solution at one point in time of one level and restricts it to the same point in time on the next coarser level.

Parameters

u – Approximate solution at a time point of one level

Returns

Input vector u at the same time point on next coarser time grid

Return type

Vector

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

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

class pymgrit.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.Mgrit(problem: List[pymgrit.core.application.Application], transfer: Optional[List[pymgrit.core.grid_transfer.GridTransfer]] = None, weight_c: float = 1.0, max_iter: int = 100, tol: float = 1e-07, nested_iteration: bool = True, cf_iter: int = 1, cycle_type: str = 'V', comm_time: mpi4py.MPI.Comm = <mpi4py.MPI.Intracomm object>, comm_space: mpi4py.MPI.Comm = <mpi4py.MPI.Comm object>, logging_lvl: int = 20, output_fcn=None, output_lvl=1, t_norm=2, random_init_guess: bool = False)

Bases: object

MGRIT solver class

Implementation of MGRIT FAS algorithm for solving time-stepping problems of the form

u_i = Phi(u_{i-1}),

where Phi propagates u_{i-1} from t = t_{i-1} to t = t_i.

It is assumed that the problems have a constant dimension and the solved space-time matrix stencil is [-Phi I].

c_exchange(lvl: int) None

Point exchange on level lvl if the first point of a process is an F-point.

Typically called after a C-point update.

Parameters

lvl – MGRIT level

c_relax(lvl: int) None

C-relaxation on level lvl.

C-relaxation updates solution values at C-points by propagating the solution from the preceeding F-points.

Parameters

lvl – MGRIT level

convergence_criterion(iteration: int) None

Stopping criterion based on the 2-norm of the space-time residual.

Computes the space-time residual

r_i = Phi(u_{i-1}) - u_i, i = 1, …. nt, r_0 = 0

Parameters

iteration – MGRIT iteration number

create_coarsest_level()

Creates vectors u and g for coarsest level

create_u(lvl: int) None

Creates solution vectors for all local time points on a given MGRIT level.

Parameters

lvl – MGRIT level

create_v_g(lvl: int) None

Creates vectors v and g for all local time points on a given MGRIT level.

Parameters

lvl – MGRIT level

error_correction(lvl: int) None

Computes the error approximation on level lvl and updates the approximation on the next finer level.

Parameters

lvl – MGRIT level

f_exchange(lvl: int) None

Point exchange on level lvl if the first point of a process is a C-point.

Typically called after an F-point update.

Parameters

lvl – MGRIT level

f_relax(lvl: int) None

F-relaxation on level lvl.

F-relaxation updates solution values at F-points by propagating the solution from one C-point to all F-points up to the next C-point.

Parameters

lvl – MGRIT level

fas_residual(lvl: int) None

Injects the fine-grid approximation and its residual from level lvl to the next coarser grid.

Parameters

lvl – MGRIT level

forward_solve(lvl: int) None

Solves the problem directly on level lvl with time stepping.

Parameters

lvl – MGRIT level

get_c_point(lvl: int) pymgrit.core.application.Application

Exchanges the first/last C-point between two processes

Parameters

lvl – MGRIT level

iteration(lvl: int, cycle_type: str, iteration: int, first_f: bool) None

MGRIT iteration on level lvl

Parameters
  • lvl – MGRIT level

  • cycle_type – Cycle type

  • iteration – Number of current iteration

  • first_f – F-relaxation at the beginning

log_info(message: str) None

Writes a message to the logger. Only one process

Parameters

message – Message

nested_iteration() None

Generates an initial approximation on the finest grid by solving the problem on the coarsest grid and interpolating the approximation to the finest level.

ouput_run_information() None

Outputs information of pyMGRIT run.

setup_points_and_comm_info(lvl: int) None

Computes local grid information for level lvl. Computes which process holds the previous and next point for each lvl and process.

Parameters

lvl – MGRIT level

solve() dict

Driver function for solving the problem using MGRIT.

Performs MGRIT iterations until a stopping criterion is fulfilled or the maximum number of iterations is reached.

Returns

dictionary with residual history, setup time, and solve time

split_into(number_points: int, number_processes: int) numpy.ndarray

Split points

Parameters
  • number_points – Number of points

  • number_processes – Number of processes

Returns

split_points(length: int, size: int, rank: int) Tuple[int, int]

Splits length points evenly in size parts and computes the index of the first point and block size of the local time interval.

Parameters
  • length – Number of points

  • size – Number of processes

  • rank – Process rank

Returns

Block size and index of first point

class pymgrit.Vector

Bases: abc.ABC

Abstract vector class for user-defined vector classes that hold information of a single time point. Every user-defined vector class must inherit from this class.

Required functions:
  • __add__

  • __sub__

  • norm

  • clone_zero

  • clone_rand

  • set_values

  • get_values

abstract clone()

Initialize vector object with same values

abstract clone_rand()

Initialize vector object with random values

abstract clone_zero()

Initialize vector object with zeros

abstract get_values(*args, **kwargs)

Get vector data

abstract norm()

Norm of a vector object

abstract pack(*args, **kwargs)

Specifying communication data

abstract set_values(*args, **kwargs)

Set vector data

abstract unpack(*args, **kwargs)

Unpacking communication data

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

pymgrit.simple_setup_problem(problem: pymgrit.core.application.Application, level: int, coarsening: int) List[pymgrit.core.application.Application]

Simple setup of a time-multigrid hierarchy for a problem.

Creates a time-multigrid hierarchy using the finest problem, the number of levels, and the coarsening factor.

Parameters
  • problem – Application problem on the finest grid

  • level – Number of time-grid levels

  • coarsening – Coarsening factor to be used for all levels

Returns

List of application problems; one application problem per time-grid level