dedalus.core.solvers
Solver classes.
Module Contents
- logger
- class SolverBase(problem, ncc_cutoff=1e-06, max_ncc_terms=None, entry_cutoff=1e-12, matrix_coupling=None, matsolver=None, bc_top=None, tau_left=None, interleave_components=None, store_expanded_matrices=None)
Base class for PDE solvers.
- Parameters
problem (Problem object) – Dedalus problem.
ncc_cutoff (float, optional) – Mode amplitude cutoff for LHS NCC expansions (default: 1e-6)
max_ncc_terms (int, optional) – Maximum terms to include in LHS NCC expansions (default: None (no limit))
entry_cutoff (float, optional) – Matrix entry cutoff to avoid fill-in from cancellation errors (default: 1e-12)
matrix_coupling (tuple of bool, optional) – Matrix coupling override.
matsolver (str or Matsolver class, optional) – Matrix solver. Default taken from config.
bc_top (bool, optional) – Whether to place boundary conditions at top of matrices. Default taken from config.
tau_left (bool, optional) – Whether to place tau columns at left of matrices. Default taken from config.
interleave_components (bool, optional) – Whether to interleave components before variables. Default taken from config.
store_expanded_matrices (bool, optional) – Whether to store right-preconditioned matrices. Default taken from config.
- build_matrices(subproblems=None, matrices=None)
Build matrices for selected subproblems.
- class EigenvalueSolver(problem, **kw)
Linear eigenvalue problem solver.
- Parameters
problem (Problem object) – Dedalus eigenvalue problem.
**kw – Other options passed to ProblemBase.
- Variables
state (list of Field objects) – Problem variables containing solution (after set_state method is called).
eigenvalues (ndarray) – Vector of eigenvalues.
eigenvectors (ndarray) – Array of eigenvectors. The eigenvector corresponding to the i-th eigenvalues is eigenvectors[:, i].
eigenvalue_subproblem (Subproblem object) – The subproblem for which the EVP has een solved.
- matsolver_default = 'MATRIX_FACTORIZER'
- matrices = ['M', 'L']
- print_subproblem_ranks(subproblems=None, target=0)
Print rank of each subproblem LHS.
- solve_dense(subproblem, rebuild_matrices=False, left=False, normalize_left=True, **kw)
Perform dense eigenvector search for selected subproblem. This routine finds all eigenvectors but is computationally expensive.
- Parameters
subproblem (Subproblem object) – Subproblem for which to solve the EVP.
rebuild_matrices (bool, optional) – Rebuild LHS matrices if coefficients have changed (default: False).
left (bool, optional) – Solve for the left eigenvectors of the system in addition to the right eigenvectors. The left eigenvectors are the right eigenvectors of the conjugate-transposed problem. Follows same definition described in scipy.linalg.eig documentation (default: False).
normalize_left (bool, optional) – 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).
**kw – Other keyword options passed to scipy.linalg.eig.
- solve_sparse(subproblem, N, target, rebuild_matrices=False, left=False, normalize_left=True, raise_on_mismatch=True, **kw)
Perform targeted sparse eigenvector search for selected subproblem. This routine finds a subset of eigenvectors near the specified target.
- Parameters
subproblem (Subproblem object) – Subproblem for which to solve the EVP.
N (int) – Number of eigenvectors to solve for. Note: the dense method may be more efficient for finding large numbers of eigenvectors.
target (complex) – Target eigenvalue for search.
rebuild_matrices (bool, optional) – Rebuild LHS matrices if coefficients have changed (default: False).
left (bool, optional) – Solve for the left eigenvectors of the system in addition to the right eigenvectors. The left eigenvectors are the right eigenvectors of the conjugate-transposed problem. Follows same definition described in scipy.linalg.eig documentation (default: False).
normalize_left (bool, optional) – 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) – Raise a RuntimeError if the left and right eigenvalues do not match (default: True).
**kw – Other keyword options passed to scipy.sparse.linalg.eig.
- set_state(index, subsystem)
Set state vector to the specified eigenmode.
- Parameters
index (int) – Index of desired eigenmode.
subsystem (Subsystem object or int) – Subsystem that will be set to the corresponding eigenmode. If an integer, the corresponding subsystem of the last specified eigenvalue_subproblem will be used.
- class LinearBoundaryValueSolver(problem, **kw)
Linear boundary value problem solver.
- Parameters
problem (Problem object) – Dedalus problem.
**kw – Other options passed to ProblemBase.
- Variables
state (list of Field objects) – Problem variables containing solution (after solve method is called).
- matsolver_default = 'MATRIX_FACTORIZER'
- matrices = ['L']
- print_subproblem_ranks(subproblems=None)
Print rank of each subproblem LHS.
- solve(subproblems=None, rebuild_matrices=False)
Solve BVP over selected subproblems.
- Parameters
subproblems (Subproblem object or list of Subproblem objects, optional) – Subproblems for which to solve the BVP (default: None (all)).
rebuild_matrices (bool, optional) – Rebuild LHS matrices if coefficients have changed (default: False).
- evaluate_handlers(handlers=None)
Evaluate specified list of handlers (all by default).
- class NonlinearBoundaryValueSolver(problem, **kw)
Nonlinear boundary value problem solver.
- Parameters
problem (Problem object) – Dedalus problem.
**kw – Other options passed to ProblemBase.
- Variables
state (list of Field objects) – Problem variables containing solution.
perturbations (list of Field objects) – Perturbations to problem variables from each Newton iteration.
iteration (int) – Current iteration.
- matsolver_default = 'MATRIX_SOLVER'
- matrices = ['dH']
- newton_iteration(damping=1)
Update solution with a Newton iteration.
- evaluate_handlers(handlers=None)
Evaluate specified list of handlers (all by default).
- class InitialValueSolver(problem, timestepper, enforce_real_cadence=100, warmup_iterations=10, **kw)
Initial value problem solver.
- Parameters
problem (Problem object) – Dedalus problem.
timestepper (Timestepper class or str) – Timestepper to use in evolving initial conditions.
enforce_real_cadence (int, optional) – Iteration cadence for enforcing Hermitian symmetry on real variables (default: 100).
warmup_iterations (int, optional) – Number of warmup iterations to disregard when computing runtime statistics (default: 10).
**kw – Other options passed to ProblemBase.
- Variables
state (list of Field objects) – Problem variables containing solution.
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.
sim_time (float) – Current simulation time.
iteration (int) – Current iteration.
dt (float) – Last timestep.
- property sim_time
- property world_time
- property wall_time
Seconds ellapsed since instantiation.
- property proceed
Check that current time and iteration pass stop conditions.
- matsolver_default = 'MATRIX_FACTORIZER'
- matrices = ['M', 'L']
- load_state(path, index=- 1, allow_missing=False)
Load state from HDF5 file. Currently can only load grid space data.
- Parameters
path (str or pathlib.Path) – Path to Dedalus HDF5 savefile
index (int, optional) – Local write index (within file) to load (default: -1)
allow_missing (bool, optional) – Do not raise an error if state variables are missing from the savefile (default: False).
- Returns
write (int) – Global write number of loaded write
dt (float) – Timestep at loaded write
- enforce_hermitian_symmetry(fields)
Transform fields to grid and back.
- step(dt)
Advance system by one iteration/timestep.
- evolve(timestep_function, log_cadence=100)
Advance system until stopping criterion is reached.
- print_subproblem_ranks(subproblems=None, dt=1)
Print rank of each subproblem LHS.
- evaluate_handlers_now(dt, handlers=None)
- evaluate_handlers(handlers=None, dt=0)
Evaluate specified list of handlers (all by default).
- log_stats(format='.4g')
Log timing statistics with specified string formatting (optional).