implementations.problem_classes.HeatEquation_ND_FD_CuPy module

class heatNd_forced(nvars=512, nu=0.1, freq=2, bc='periodic', order=2, stencil_type='center', lintol=1e-12, liniter=10000, solver_type='direct')[source]

Bases: Problem

This class implements the unforced \(N\)-dimensional heat equation with periodic boundary conditions

\[\frac{\partial u}{\partial t} = \nu \left( \frac{\partial^2 u}{\partial x^2_1} + .. + \frac{\partial^2 u}{\partial x^2_N} \right)\]

for \((x_1,..,x_N) \in [0, 1]^{N}\) with \(N \leq 3\). The initial solution is of the form

\[u({\bf x},0) = \prod_{i=1}^N \sin(\pi k_i x_i).\]

The spatial term is discretized using finite differences.

This class uses the CuPy package in order to make pySDC available for GPUs.

Parameters:
  • nvars (int, optional) – Spatial resolution (same in all dimensions). Using a tuple allows to consider several dimensions, e.g nvars=(16,16) for a 2D problem.

  • nu (float, optional) – Diffusion coefficient \(\nu\).

  • freq (int, optional) – Spatial frequency \(k_i\) of the initial conditions, can be tuple.

  • stencil_type (str, optional) – Type of the finite difference stencil.

  • order (int, optional) – Order of the finite difference discretization.

  • lintol (float, optional) – Tolerance for spatial solver.

  • liniter (int, optional) – Max. iterations number for spatial solver.

  • solver_type (str, optional) – Solve the linear system directly or using CG.

  • bc (str, optional) – Boundary conditions, either 'periodic' or 'dirichlet'.

  • sigma (float, optional) –

    If freq=-1 and ndim=1, uses a Gaussian initial solution of the form

    \[u(x,0) = e^{ \frac{\displaystyle 1}{\displaystyle 2} \left( \frac{\displaystyle x-1/2}{\displaystyle \sigma} \right)^2 }\]

A

FD discretization matrix of the ND grad operator.

Type:

sparse matrix (CSC)

Id

Identity matrix of the same dimension as A

Type:

sparse matrix (CSC)

dtype_f

alias of imex_cupy_mesh

dtype_u

alias of cupy_mesh

eval_f(u, t)[source]

Routine to evaluate the right-hand side of the problem.

Parameters:
  • u (dtype_u) – Current values of the numerical solution.

  • t (float) – Current time of the numerical solution is computed.

Returns:

f – The right-hand side of the problem.

Return type:

dtype_f

property ndim

Number of dimensions of the spatial problem

solve_system(rhs, factor, u0, t)[source]

Simple linear solver for \((I-factor\cdot A)\vec{u}=\vec{rhs}\).

Parameters:
  • rhs (dtype_f) – Right-hand side for the linear system.

  • factor (float) – Abbrev. for the local stepsize (or any other factor required).

  • u0 (dtype_u) – Initial guess for the iterative solver.

  • t (float) – Current time (e.g. for time-dependent BCs).

Returns:

me – Solution.

Return type:

dtype_u

u_exact(t)[source]

Routine to compute the exact solution at time \(t\).

Parameters:

t (float) – Time of the exact solution.

Returns:

me – Exact solution.

Return type:

dtype_u

class heatNd_unforced(nvars=512, nu=0.1, freq=2, bc='periodic', order=2, stencil_type='center', lintol=1e-12, liniter=10000, solver_type='direct')[source]

Bases: heatNd_forced

This class implements the forced \(N\)-dimensional heat equation with periodic boundary conditions

\[\frac{\partial u}{\partial t} = \nu \left( \frac{\partial^2 u}{\partial x^2_1} + .. + \frac{\partial^2 u}{\partial x^2_N} \right) + f({\bf x}, t)\]

for \((x_1,..,x_N) \in [0, 1]^{N}\) with \(N \leq 3\), and forcing term

\[f({\bf x}, t) = \prod_{i=1}^N \sin(\pi k_i x_i) \left( \nu \pi^2 \sum_{i=1}^N k_i^2 \cos(t) - \sin(t) \right),\]

where \(k_i\) denotes the frequency in the \(i^{th}\) dimension. The exact solution is

\[u({\bf x}, t) = \prod_{i=1}^N \sin(\pi k_i x_i) \cos(t).\]

The spatial term is discretized using finite differences.

The implementation is this class uses the CuPy package in order to make pySDC available for GPUs.

dtype_f

alias of cupy_mesh

eval_f(u, t)[source]

Routine to evaluate the right-hand side of the problem.

Parameters:
  • u (dtype_u) – Current values of the numerical solution.

  • t (float) – Current time of the numerical solution is computed.

Returns:

f – The right-hand side of the problem.

Return type:

dtype_f

u_exact(t)[source]

Routine to compute the exact solution at time \(t\).

Parameters:

t (float) – Time of the exact solution.

Returns:

me – The exact solution.

Return type:

dtype_u