core.Problem module¶
Description¶
Module containing the base Problem class for pySDC
- class WorkCounter[source]¶
Bases:
object
Utility class for counting iterations.
Contains one attribute niter initialized to zero during instantiation, which can be incremented by calling object as a function, e.g
>>> count = WorkCounter() # => niter = 0 >>> count() # => niter = 1 >>> count() # => niter = 2
- class ptype(init)[source]¶
Bases:
RegisterParams
Prototype class for problems, just defines the attributes essential to get started.
- Parameters:
init (list of args) – Argument(s) used to initialize data types.
dtype_u (type) – Variable data type. Should generate a data variable using dtype_u(init).
dtype_f (type) – RHS data type. Should generate a data variable using dtype_f(init).
- logger¶
custom logger for problem-related logging.
- Type:
logging.Logger
- dtype_f = None¶
- dtype_u = None¶
- eval_f(u, t)[source]¶
Abstract interface to RHS computation of the ODE
- Parameters:
u (dtype_u) – Current values.
t (float) – Current time.
- Returns:
f – The RHS values.
- Return type:
dtype_f
- property f_init¶
Generate a data variable for RHS
- generate_scipy_reference_solution(eval_rhs, t, u_init=None, t_init=None, **kwargs)[source]¶
Compute a reference solution using scipy.solve_ivp with very small tolerances. Keep in mind that scipy needs the solution to be a one dimensional array. If you are solving something higher dimensional, you need to make sure the function eval_rhs takes a flattened one-dimensional version as an input and output, but reshapes to whatever the problem needs for evaluation.
The keyword arguments will be passed to scipy.solve_ivp. You should consider passing method=’BDF’ for stiff problems and to accelerate that you can pass a function that evaluates the Jacobian with arguments jac(t, u) as jac=jac.
- Parameters:
eval_rhs (function) – Function evaluate the full right hand side. Must have signature eval_rhs(float: t, numpy.1darray: u)
t (float) – current time
u_init (pySDC.implementations.problem_classes.Lorenz.dtype_u) – initial conditions for getting the exact solution
t_init (float) – the starting time
- Returns:
Reference solution
- Return type:
numpy.ndarray
- get_fig()[source]¶
Get a figure suitable to plot the solution of this problem
- Returns:
self.fig
- Return type:
matplotlib.pyplot.figure.Figure
- logger = <Logger problem (WARNING)>¶
- plot(u, t=None, fig=None)[source]¶
Plot the solution. Please supply a figure with the same structure as returned by
self.get_fig
.- Parameters:
u (dtype_u) – Solution to be plotted
t (float) – Time to display at the top of the figure
fig (matplotlib.pyplot.figure.Figure) – Figure with the correct structure
- Return type:
None
- property u_init¶
Generate a data variable for u