dedalus.core.solvers

Classes for solving differential equations.

Module Contents

logger
class EigenvalueSolver(problem, matsolver=None)

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(pencil, rebuild_coeffs=False, normalize_left=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)

  • normalize_left (bool, optional) – Flag to normalize the left eigenvectors (if left=True) such that the modified left eigenvectors form a biorthonormal (not just biorthogonal) set with respect to the right eigenvectors. (default: True)

  • Other keyword options passed to scipy.linalg.eig.

solve_sparse(pencil, N, target, rebuild_coeffs=False, left=False, normalize_left=True, raise_on_mismatch=True, **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)

  • left (bool, optional) – Flag to solve for the left eigenvectors of the system (IN ADDITION TO the right eigenvectors), defined as the right eigenvectors of the conjugate-transposed problem. Follows same definition described in scipy.linalg.eig documentation. (default: False)

  • normalize_left (bool, optional) – Flag to normalize the left eigenvectors such that the modified left eigenvectors form a biorthonormal (not just biorthogonal) set with respect to the right eigenvectors. (default: True)

  • raise_on_mismatch (bool, optional) – Flag to raise a RuntimeError if the eigenvalues of the conjugate- transposed problem (i.e. the left eigenvalues) do not match the original (or “right”) eigenvalues. (default: True)

  • Other keyword options passed to scipy.sparse.linalg.eigs.

set_state(index, left=False, modified_left=False)

Set state vector to the specified eigenmode.

Parameters
  • index (int) – Index of desired eigenmode

  • left (bool, optional) – If true, sets state vector to a left (or adjoint) eigenmode instead of right eigenmode unless modified_left=True (default: False)

  • modified_left (bool, optional) – If true, sets state vector to a modified left eigenmode, which is dual (under the standard complex dot product in coefficient space) to the corresponding right eigenmode. Supersedes left=True (default: False)

solve(*args, **kw)

Deprecated. Use solve_dense instead.

class LinearBoundaryValueSolver(problem, matsolver=None)

Linear boundary value problem solver.

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

  • matsolver (matsolver class or name, optional) – Matrix solver routine (default set by config file).

Variables

state (system object) – System containing solution fields (after solve method is called).

solve(rebuild_coeffs=False)

Solve BVP.

class NonlinearBoundaryValueSolver(problem, matsolver=None)

Nonlinear boundary value problem solver.

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

  • matsolver (matsolver class or name, optional) – Matrix solver routine (default set by config file).

Variables

state (system object) – System containing solution fields (after solve method is called)

newton_iteration(damping=1)

Update solution with a Newton iteration.

class InitialValueSolver(problem, timestepper, matsolver=None, enforce_real_cadence=100, warmup_iterations=10)

Initial value problem solver.

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

  • timestepper (timestepper class or name) – Timestepper to use in evolving initial conditions

  • matsolver (matsolver class or name, optional) – Matrix solver routine (default set by config file).

  • enforce_real_cadence (int, optional) – Iteration cadence for enforcing Hermitian symmetry on real variables (default: 100).

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

property sim_time
property proceed

Check that current time and iteration pass stop conditions.

property ok

Deprecated. Use ‘solver.proceed’.

get_world_time()
load_state(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()

Build array of finite handler sim_dt cadences.

step(dt, trim=False)

Advance system by one iteration/timestep.

evolve(timestep_function)

Advance system until stopping criterion is reached.

evaluate_handlers_now(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.

log_stats(format='.4g')

Log timing statistics with specified string formatting (optional).