implementations.convergence_controller_classes.estimate_embedded_error module

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

Bases: ConvergenceController

The embedded error is obtained by computing two solutions of different accuracy and pretending the more accurate one is an exact solution from the point of view of the less accurate solution. In practice, we like to compute the solutions with different order methods, meaning that in SDC we can just subtract two consecutive sweeps, as long as you make sure your preconditioner is compatible, which you have to just try out…

dependencies(controller, description, **kwargs)[source]

Load the convergence controller that stores the solution of the last sweep unless we are doing Runge-Kutta. Add the hook for recording the error.

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

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

Returns:

None

estimate_embedded_error_serial(L)[source]

Estimate the serial embedded error, which may need to be modified for a parallel estimate.

Depending on the type of sweeper, the lower order solution is stored in a different place.

Parameters:

L (pySDC.level) – The level

Returns:

The embedded error estimate

Return type:

dtype_u

classmethod get_implementation(flavor='standard', useMPI=False)[source]

Retrieve the implementation for a specific flavor of this class.

Parameters:

flavor (str) – The implementation that you want

Returns:

The child class that implements the desired flavor

Return type:

cls

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

Estimate the local error here.

If you are doing MSSDC, this is the global error within the block in Gauss-Seidel mode. In Jacobi mode, I haven’t thought about what this is.

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

  • S (pySDC.Step) – The current step

Returns:

None

reset_status_variables(controller, **kwargs)[source]

Reset status variables. This is called in the restart_block function. :param controller: The controller :type controller: pySDC.Controller

Returns:

None

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

Add a default value for control order to the parameters and check if we are using a Runge-Kutta sweeper

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

  • params (dict) – Parameters for the convergence controller

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

Returns:

Updated parameters

Return type:

dict

setup_status_variables(controller, **kwargs)[source]

Add the embedded error variable to the error function.

Parameters:

controller (pySDC.Controller) – The controller

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

Bases: ConvergenceController

Estimates an embedded error based on changing the underlying quadrature rule. The error estimate is stored as error_embedded_estimate_collocation in the status of the level. Note that we only compute the estimate on the finest level. The error is stored as a tuple with the first index denoting to which iteration it belongs. This is useful since the error estimate is not available immediately after, but only when the next collocation problem is converged to make sure the two solutions are of different accuracy.

Changing the collocation method between iterations happens using the AdaptiveCollocation convergence controller. Please refer to that for documentation on how to use this. Just pass the parameters for that convergence controller as adaptive_coll_params to the parameters for this one and they will be passed on when the AdaptiveCollocation convergence controller is automatically added while loading dependencies.

dependencies(controller, description, **kwargs)[source]

Load the AdaptiveCollocation convergence controller to switch between collocation problems between iterations.

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

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

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

Compute the embedded error as the difference between the interpolated and the current solution on the finest level.

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

  • step (pySDC.Step.step) – The current step

reset_status_variables(controller, **kwargs)[source]

Reset status variables. This is called in the restart_block function. :param controller: The controller :type controller: pySDC.Controller

Returns:

None

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

Add a default value for control order to the parameters

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

  • params (dict) – Parameters for the convergence controller

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

Returns:

Updated parameters

Return type:

dict

setup_status_variables(controller, **kwargs)[source]

Add the embedded error variable to the levels and add a status variable for previous steps.

Parameters:

controller (pySDC.Controller) – The controller

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

Bases: EstimateEmbeddedError

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

Compute embedded error estimate on the last node of each level In serial this is the local error, but in block Gauss-Seidel MSSDC this is a semi-global error in each block

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

  • S (pySDC.Step) – The current step

Returns:

None

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

Bases: EstimateEmbeddedError

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

Compute embedded error estimate on the last node of each level In serial this is the local error, but in block Gauss-Seidel MSSDC this is a semi-global error in each block

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

  • S (pySDC.Step) – The current step

Returns:

None

reset_buffers_nonMPI(controller, **kwargs)[source]

Reset buffers for imitated communication.

Parameters:

controller (pySDC.controller) – The controller

Returns:

None

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

Add option for averaging the local errors.

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

  • params (dict) – Parameters for the convergence controller

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

Returns:

Updated parameters

Return type:

dict