dedalus.core.solvers

Classes for solving differential equations.

Module Contents

PERMC_SPEC
USE_UMFPACK
logger
class EigenvalueSolver(problem)

Solves linear eigenvalue problems for oscillation frequency omega, (d_t -> -i omega) for a given pencil, and stores the eigenvalues and eigenvectors. The set_state method can be used to set solver.state to the specified eigenmode.

Parameters:

problem (problem object) – Problem describing system of differential equations and constraints

Variables:
  • state (system object) – System containing solution fields (after solve method is called)
  • eigenvalues (numpy array) – Contains a list of eigenvalues omega
  • eigenvectors (numpy array) – Contains a list of eigenvectors. The eigenvector corresponding to the ith eigenvalue is in eigenvectors[:,i]
  • eigenvalue_pencil (pencil) – The pencil for which the eigenvalue problem has been solved.
solve_dense(self, pencil, rebuild_coeffs=False, **kw)

Solve EVP for selected pencil.

Parameters:
  • pencil (pencil object) – Pencil for which to solve the EVP
  • rebuild_coeffs (bool, optional) – Flag to rebuild cached coefficient matrices (default: False)
  • Other keyword options passed to scipy.linalg.eig.
solve_sparse(self, pencil, N, target, rebuild_coeffs=False, **kw)

Perform targeted sparse eigenvalue search for selected pencil.

Parameters:
  • pencil (pencil object) – Pencil for which to solve the EVP
  • N (int) – Number of eigenmodes to solver for. Note: the dense method may be more efficient for finding large numbers of eigenmodes.
  • target (complex) – Target eigenvalue for search.
  • rebuild_coeffs (bool, optional) – Flag to rebuild cached coefficient matrices (default: False)
  • Other keyword options passed to scipy.sparse.linalg.eigs.
set_state(self, index)

Set state vector to the specified eigenmode.

Parameters:index (int) – Index of desired eigenmode
solve(self, *args, **kw)

Deprecated. Use solve_dense instead.

class LinearBoundaryValueSolver(problem)

Linear boundary value problem solver.

Parameters:problem (problem object) – Problem describing system of differential equations and constraints
Variables:state (system object) – System containing solution fields (after solve method is called)
solve(self, rebuild_coeffs=False)

Solve BVP.

class NonlinearBoundaryValueSolver(problem)

Nonlinear boundary value problem solver.

Parameters:problem (problem object) – Problem describing system of differential equations and constraints
Variables:state (system object) – System containing solution fields (after solve method is called)
newton_iteration(self, damping=1)

Update solution with a Newton iteration.

class InitialValueSolver(problem, timestepper)

Initial value problem solver.

Parameters:
  • problem (problem object) – Problem describing system of differential equations and constraints
  • timestepper (timestepper class) – Timestepper to use in evolving initial conditions
Variables:
  • state (system object) – System containing current solution fields
  • dt (float) – Timestep
  • stop_sim_time (float) – Simulation stop time, in simulation units
  • stop_wall_time (float) – Wall stop time, in seconds from instantiation
  • stop_iteration (int) – Stop iteration
  • time (float) – Current simulation time
  • iteration (int) – Current iteration
sim_time
ok

Check that current time and iteration pass stop conditions.

get_world_time(self)
load_state(self, path, index=-1)

Load state from HDF5 file.

Parameters:
  • path (str or pathlib.Path) – Path to Dedalus HDF5 savefile
  • index (int, optional) – Local write index (within file) to load (default: -1)
Returns:

  • write (int) – Global write number of loaded write
  • dt (float) – Timestep at loaded write

sim_dt_cadences(self)

Build array of finite handler sim_dt cadences.

step(self, dt, trim=False)

Advance system by one iteration/timestep.

evolve(self, timestep_function)

Advance system until stopping criterion is reached.

evaluate_handlers_now(self, dt, handlers=None)

Evaluate all handlers right now. Useful for writing final outputs.

by default, all handlers are evaluated; if a list is given only those will be evaluated.