4.3.2 pde.pdes.base module
Base class for defining partial differential equations.
- class PDEBase(*, noise=0, rng=None)[source]
Bases:
object
Base class for defining partial differential equations (PDEs)
Custom PDEs can be implemented by subclassing
PDEBase
to specify the evolution rate. In the simple case of deterministic PDEs, the methodsPDEBase.evolution_rate()
andPDEBase._make_pde_rhs_numba()
need to be overwritten for supporting the numpy and numba backend, respectively.- Parameters:
noise (float or
ndarray
) – Variance of the additive Gaussian white noise that is supported for all PDEs by default. If set to zero, a deterministic partial differential equation will be solved. Different noise magnitudes can be supplied for each field in coupled PDEs.rng (
Generator
) – Random number generator (default:default_rng()
) used for stochastic simulations. Note that this random number generator is only used for numpy function, while compiled numba code uses the random number generator of numba. Moreover, in simulations using multiprocessing, setting the same generator in all processes might yield unintended correlations in the simulation results.
Note
If more complicated noise structures are required, the methods
PDEBase.noise_realization()
andPDEBase._make_noise_realization_numba()
need to be overwritten for the numpy and numba backend, respectively.- cache_rhs: bool = False
Flag indicating whether the right hand side of the equation should be cached. If True, the same implementation is used in subsequent calls to solve. Note that the cache is only invalidated when the grid of the underlying state changes. Consequently, the simulation might lead to wrong results if the parameters of the PDE are changed after the first call. This option is thus disabled by default and should be used with care.
- Type:
- check_implementation: bool = True
Flag determining whether numba-compiled functions should be checked against their numpy counter-parts. This can help with implementing a correct compiled version for a PDE class.
- Type:
- check_rhs_consistency(state, t=0, *, tol=1e-07, rhs_numba=None, **kwargs)[source]
Check the numba compiled right hand side versus the numpy variant.
- Parameters:
state (
FieldBase
) – The state for which the evolution rates should be comparedt (float) – The associated time point
tol (float) – Acceptance tolerance. The check passes if the evolution rates differ by less then this value
rhs_numba (callable) – The implementation of the numba variant that is to be checked. If omitted, an implementation is obtained by calling
PDEBase._make_pde_rhs_numba_cached()
.
- Return type:
None
- complex_valued: bool = False
Flag indicating whether the right hand side is a complex-valued PDE, which requires all involved variables to have complex data type.
- Type:
- explicit_time_dependence: bool | None = None
Flag indicating whether the right hand side of the PDE has an explicit time dependence.
- Type:
- property is_sde: bool
flag indicating whether this is a stochastic differential equation
The
BasePDF
class supports additive Gaussian white noise, whose magnitude is controlled by the noise property. In this case, is_sde is True if self.noise != 0.- Type:
- make_pde_rhs(state, backend='auto', **kwargs)[source]
Return a function for evaluating the right hand side of the PDE.
- Parameters:
state (
FieldBase
) – An example for the state from which the grid and other information can be extracted.backend (str) – Determines how the function is created. Accepted values are ‘numpy’ and ‘numba’. Alternatively, ‘auto’ lets the code pick the optimal backend.
- Returns:
Function determining the right hand side of the PDE
- Return type:
callable
- make_post_step_hook(state)[source]
Returns a function that is called after each step.
This function receives three arguments: the current state as a numpy array, the current time point, and a numpy array that can store data for the hook function. The function can modify the state data in place. If the function makes use of the data feature, it must replace the data in place.
The hook can also be used to abort the simulation when a user-defined condition is met by raising StopIteration. Note that this interrupts the inner-most loop, so that some final information might be still reflect the values they assumed at the last tracker interrupt. Additional information (beside the current state) should be returned by the post_step_data.
Example
The following code provides an example that creates a hook function that limits the state to a maximal value of 1 and keeps track of the total correction that is applied. This is achieved using post_step_data, which is initialized with the second value (0) returned by the method and incremented each time the hook is called.
def make_post_step_hook(self, state): def post_step_hook(state_data, t, post_step_data): i = state_data > 1 # get violating entries overshoot = (state_data[i] - 1).sum() # get total correction state_data[i] = 1 # limit data entries post_step_data += overshoot # accumulate total correction return post_step_hook, 0.0 # hook function and initial value
- Parameters:
state (
FieldBase
) – An example for the state from which the grid and other information can be extracted- Returns:
- The first entry is the function that implements the hook. The second
entry gives the initial data that is used as auxiliary data in the hook. This can be None if no data is used.
- Return type:
- make_sde_rhs(state, backend='auto', **kwargs)[source]
Return a function for evaluating the right hand side of the SDE.
- Parameters:
state (
FieldBase
) – An example for the state from which information can be extractedbackend (str) – Determines how the function is created. Accepted values are ‘numpy’ and ‘numba’. Alternatively, ‘auto’ lets the code pick the optimal backend.
- Returns:
Function determining the deterministic part of the right hand side of the PDE together with a noise realization.
- Return type:
- noise_realization(state, t=0, *, label='Noise realization')[source]
Returns a realization for the noise.
- solve(state, t_range, dt=None, tracker='auto', *, solver='explicit', ret_info=False, **kwargs)[source]
Solves the partial differential equation.
The method constructs a suitable solver (
SolverBase
) and controller (Controller
) to advance the state over the temporal range specified by t_range. This method only exposes the most common functions, so explicit construction of these classes might offer more flexibility.- Parameters:
state (
FieldBase
) – The initial state (which also defines the spatial grid).t_range (float or tuple) – Sets the time range for which the PDE is solved. This should typically be a tuple of two numbers, (t_start, t_end), specifying the initial and final time of the simulation. If only a single value is given, it is interpreted as t_end and the time range is (0, t_end).
dt (float) – Time step of the chosen stepping scheme. If None, a default value based on the stepper will be chosen. If an adaptive stepper is used (supported by
ScipySolver
andExplicitSolver
), dt sets the initial time step.tracker (TrackerCollectionDataType) – Defines trackers that process the state of the simulation at specified times. A tracker is either an instance of
TrackerBase
or a string identifying a tracker (possible identifiers can be obtained by callingget_named_trackers()
). Multiple trackers can be specified as a list. The default value auto checks the state for consistency (tracker ‘consistency’) and displays a progress bar (tracker ‘progress’) whentqdm
is installed. More general trackers are defined intrackers
, where all options are explained in detail. In particular, the time points where the tracker analyzes data can be chosen when creating a tracker object explicitly.solver (
SolverBase
or str) – Specifies the method for solving the differential equation. This can either be an instance ofSolverBase
or a descriptive name like ‘explicit’ or ‘scipy’. The valid names are given bypde.solvers.registered_solvers()
. Details of the solvers and additional features (like adaptive time steps) are explained insolvers
.ret_info (bool) – Flag determining whether diagnostic information about the solver process should be returned. Note that the same information is also available as the
diagnostics
attribute.**kwargs – Additional keyword arguments are forwarded to the solver class chosen with the solver argument. In particular,
ExplicitSolver
supports several schemes and an adaptive stepper can be enabled usingadaptive=True
. Conversely,ScipySolver
accepts the additional arguments ofscipy.integrate.solve_ivp()
.
- Returns:
The state at the final time point. If ret_info == True, a tuple with the final state and a dictionary with additional information is returned. Note that None instead of a field is returned in multiprocessing simulations if the current node is not the main MPI node.
- Return type: