implementations.convergence_controller_classes.estimate_extrapolation_error module

class EstimateExtrapolationErrorBase(controller, params, description, **kwargs)[source]

Bases: ConvergenceController

Abstract base class for extrapolated error estimates

This error estimate extrapolates a solution based on Taylor expansions using solutions of previous time steps. In particular, child classes need to implement how to make these solutions available, which works differently for MPI and non-MPI versions.

check_parameters(controller, params, description, **kwargs)[source]

Check whether parameters are compatible with whatever assumptions went into the step size functions etc.

Parameters:
  • controller (pySDC.Controller) – The controller

  • params (dict) – The params passed for this specific convergence controller

  • description (dict) – The description object used to instantiate the controller

Returns:

Whether the parameters are compatible str: Error message

Return type:

bool

get_extrapolation_coefficients(t, dt, t_eval)[source]

This function solves a linear system where in the matrix A, the row index reflects the order of the derivative in the Taylor expansion and the column index reflects the particular step and whether its u or f from that step. The vector b on the other hand, contains a 1 in the first entry and zeros elsewhere, since we want to compute the value itself and all the derivatives should vanish after combining the Taylor expansions. This works to the order of the number of rows and since we want a square matrix for solving, we need the same amount of columns, which determines the memory overhead, since it is equal to the solutions / rhs that we need in memory at the time of evaluation.

This is enough to get the extrapolated solution, but if we want to compute the local error, we have to compute a prefactor. This is based on error accumulation between steps (first step’s solution is exact plus 1 LTE, second solution is exact plus 2 LTE and so on), which can be computed for adaptive step sizes as well. However, this is only true for linear problems, which means we expect the error estimate to work less well for non-linear problems.

Since only time differences are important for computing the coefficients, we need to compute this only once when using constant step sizes. When we allow the step size to change, however, we need to recompute this in every step, which is activated by the use_adaptivity parameter.

Solving for the coefficients requires solving a dense linear system of equations. The number of unknowns is equal to the order of the Taylor expansion, so this step should be cheap compared to the solves in each SDC iteration.

The function stores the computed coefficients in the self.coeff variables.

Parameters:
  • t (list) – The list of times at which we have solutions available

  • dt (list) – The step sizes used for computing these solutions (needed for the prefactor)

  • t_eval (float) – The time we want to extrapolate to

Returns:

None

reset_status_variables(controller, **kwargs)[source]

Add variable for extrapolated error

Parameters:

controller (pySDC.Controller) – The controller

Returns:

None

setup(controller, params, description, **kwargs)[source]

The extrapolation based method requires storage of previous values of u, f, t and dt and also requires solving a linear system of equations to compute the Taylor expansion finite difference style. Here, all variables are initialized which are needed for this process.

Parameters:
  • controller (pySDC.Controller) – The controller

  • params (dict) – The params passed for this specific convergence controller

  • description (dict) – The description object used to instantiate the controller

Returns:

Updated parameters with default values

Return type:

dict

setup_status_variables(controller, **kwargs)[source]

Initialize coefficient variables and add variable to the levels for extrapolated error

Parameters:

controller (pySDC.controller) – The controller

Returns:

None

store_values(S, **kwargs)[source]

Store the required attributes of the step to do the extrapolation. We only care about the last collocation node on the finest level at the moment.

Parameters:

S (pySDC.Step) – The current step

Returns:

None

class EstimateExtrapolationErrorNonMPI(controller, params, description, **kwargs)[source]

Bases: EstimateExtrapolationErrorBase

Implementation of the extrapolation error estimate for the non-MPI controller.

get_extrapolated_error(S, **kwargs)[source]

The extrapolation estimate combines values of u and f from multiple steps to extrapolate and compare to the solution obtained by the time marching scheme.

Parameters:

S (pySDC.Step) – The current step

Returns:

None

get_extrapolated_solution(S, **kwargs)[source]

Combine values from previous steps to extrapolate.

Parameters:

S (pySDC.Step) – The current step

Returns:

The extrapolated solution

Return type:

dtype_u

post_iteration_processing(controller, S, **kwargs)[source]
We perform three key operations here in the last iteration:
  • Compute the error estimate

  • Compute the coefficients if needed

  • Store the values of the step if we pretend not to for the no-memory overhead version

Parameters:
  • controller (pySDC.Controller) – The controller

  • S (pySDC.Step) – The current step

Returns:

None

prepare_next_block(controller, S, size, time, Tend, MS, **kwargs)[source]

If the no-memory-overhead version is used, we need to delete stuff that shouldn’t be available. Otherwise, we need to store all the stuff that we can.

Parameters:
  • controller (pySDC.Controller) – The controller

  • S (pySDC.step) – The current step

  • size (int) – Number of ranks

  • time (float) – The current time

  • Tend (float) – The final time

  • MS (list) – Active steps

Returns:

None

setup(controller, params, description, **kwargs)[source]

Add a no parameter ‘no_storage’ which decides whether the standard or the no-memory-overhead version is run, where only values are used for extrapolation which are in memory of other processes

Parameters:
  • controller (pySDC.Controller) – The controller

  • params (dict) – The params passed for this specific convergence controller

  • description (dict) – The description object used to instantiate the controller

Returns:

Updated parameters with default values

Return type:

dict

setup_status_variables(controller, **kwargs)[source]

Initialize storage variables.

Parameters:

controller (pySDC.controller) – The controller

Returns:

None

class EstimateExtrapolationErrorWithinQ(controller, params, description, **kwargs)[source]

Bases: EstimateExtrapolationErrorBase

This convergence controller estimates the local error based on comparing the SDC solution to an extrapolated solution within the quadrature matrix. Collocation methods compute a high order solution from a linear combination of solutions at intermediate time points. While the intermediate solutions (a.k.a. stages) don’t share the order of accuracy with the solution at the end of the interval, for SDC we know that the order is equal to the number of nodes + 1 (locally). This is because the solution to the collocation problem is a polynomial approximation of order of the number of nodes. That means we can do a Taylor expansion around the end point of the interval to higher order and after cancelling terms just like we are used to with the extrapolation based error estimate across multiple steps, we get an error estimate that is of the order of accuracy of the stages. This can be used for adaptivity, for instance, with the nice property that it doesn’t matter how we arrived at the converged collocation solution, as long as we did. We don’t rely on knowing the order of accuracy after every sweep, only after convergence of the collocation problem has been achieved, which we can check from the residual.

post_iteration_processing(controller, S, **kwargs)[source]

Compute the extrapolated error estimate here if the step is converged.

Parameters:
  • controller (pySDC.Controller) – The controller

  • S (pySDC.Step) – The current step

Returns:

None

setup(controller, params, description, **kwargs)[source]

We need this convergence controller to become active after the check for convergence, because we need the step to be converged.

Parameters:
  • controller (pySDC.Controller) – The controller

  • params (dict) – The params passed for this specific convergence controller

  • description (dict) – The description object used to instantiate the controller

Returns:

Updated parameters with default values

Return type:

dict