eigentools.eigenproblem

Module Contents

logger
class Eigenproblem(EVP, reject=True, factor=1.5, scales=1, drift_threshold=1000000.0, use_ordinal=False, grow_func=lambda x: ..., freq_func=lambda x: ...)

An object for feature-rich eigenvalue analysis.

Eigenproblem provides support for common tasks in eigenvalue analysis. Dedalus EVP objects compute raw eigenvalues and eigenvectors for a given problem; Eigenproblem provides support for numerous common tasks required for scientific use of those solutions. This includes rejection of inaccurate eigenvalues and analysis of those rejection criteria, plotting of eigenmodes and spectra, and projection of 1-D eigenvectors onto 2- or 3-D domains for use as initial conditions in subsequent initial value problems.

Additionally, Eigenproblems can compute epsilon-pseudospectra for arbitrary Dedalus differential-algebraic equations.

Parameters
  • EVP (dedalus.core.problems.EigenvalueProblem) – The Dedalus EVP object containing the equations to be solved

  • reject (bool, optional) – whether or not to reject spurious eigenvalues (default: True)

  • factor (float, optional) – The factor by which to multiply the resolution. NB: this must be a rational number such that factor times the resolution of EVP is an integer. (default: 1.5)

  • scales (float, optional) – A multiple for setting the grid resolution. (default: 1)

  • drift_threshold (float, optional) – Inverse drift ratio threshold for keeping eigenvalues during rejection (default: 1e6)

  • use_ordinal (bool, optional) – If true, use ordinal method from Boyd (1989); otherwise use nearest (default: False)

  • grow_func (func) – A function that takes a complex input and returns the growth rate as defined by the EVP (default: uses real part)

  • freq_func (func) – A function that takes a complex input and returns the frequency as defined by the EVP (default: uses imaginary part)

Variables
  • evalues (ndarray) – Lists “good” eigenvalues

  • evalues_low (ndarray) – Lists eigenvalues from low resolution solver (i.e. the resolution of the specified EVP)

  • evalues_high (ndarray) – Lists eigenvalues from high resolution solver (i.e. factor times specified EVP resolution)

  • pseudospectrum (ndarray) – epsilon-pseudospectrum computed at specified points in the complex plane

  • ps_real (ndarray) – real coordinates for epsilon-pseudospectrum

  • ps_imag (ndarray) – imaginary coordinates for epsilon-pseudospectrum

Notes

See references for algorithms in individual method docstrings.

grid(self)

get grid points for eigenvectors.

solve(self, sparse=False, parameters=None, pencil=0, N=15, target=0, **kwargs)

solve underlying eigenvalue problem.

Parameters
  • sparse (bool, optional) – If true, use sparse solver, otherwise use dense solver (default: False)

  • parameters (dict, optional) – A dict giving parameter names and values to the EVP. If None, use values specified at EVP construction time. (default: None)

  • pencil (int, optional) – The EVP pencil to be solved. (default: 0)

  • N (int, optional) – The number of eigenvalues to find if using a sparse solver (default: 15)

  • target (complex, optional) – The target value to search for when using sparse solver (default: 0+0j)

eigenmode(self, index, scales=None, all_modes=False)

Returns Dedalus FieldSystem object containing the eigenmode given by index.

Parameters
  • index (int) – index of eigenvalue corresponding to desired eigenvector

  • scales (float) – A multiple for setting the grid resolution. If not None, will overwrite self.scales. (default: None)

  • all_modes (bool, optional) – If True, index specifies the unsorted index of the low-resolution EVP; otherwise it is the index corresponding to the self.evalues order (default: False)

growth_rate(self, parameters=None, **kwargs)

returns the maximum growth rate, defined by self.grow_func(), the index of the maximal mode, and the frequency of that mode. If there is no growing mode, returns the slowest decay rate.

also returns the index of the fastest growing mode. If there are no good eigenvalues, returns np.nan for all three quantities.

Returns

growth_rate, index, freqency

Return type

tuple of ints

plot_mode(self, index, fig_height=8, norm_var=None, scales=None, all_modes=False)

plots eigenvector corresponding to specified index.

By default, the plot will show the real and complex parts of the unnormalized components of the eigenmode. If a norm_var is specified, all components will be scaled such that variable chosen is purely real and has unit amplitude.

Parameters
  • index (int) – index of eigenvalue corresponding to desired eigenvector

  • fig_height (float, optional) – Height of constructed figure (default: 8)

  • norm_var (str) – If not None, selects the field in the eigenmode with which to normalize. Otherwise, plots the unnormalized eigenmode. (default: None)

  • scales (float) – A multiple for setting the grid resolution. If not None, will overwrite self.scales. (default: None)

  • all_modes (bool, optional) – If True, index specifies the unsorted index of the low-resolution EVP; otherwise it is the index corresponding to the self.evalues order (default: False)

Returns

Return type

matplotlib.figure.Figure

project_mode(self, index, domain, transverse_modes, all_modes=False)

projects a mode specified by index onto a domain of higher dimension.

Parameters
  • index – an integer giving the eigenmode to project

  • domain – a domain to project onto

  • transverse_modes – a tuple of mode numbers for the transverse directions

Returns

Return type

dedalus.core.system.FieldSystem

write_global_domain(self, field_system, base_name='IVP_output')

Given a field system, writes a Dedalus HDF5 file.

Typically, one would use this to write a field system constructed by project_mode.

Parameters
  • field_system (dedalus.core.system.FieldSystem) – A field system containing the data to be written

  • base_name (str, optional) – The base filename of the resulting HDF5 file. (default: IVP_output)

calc_ps(self, k, zgrid, mu=0.0, pencil=0, inner_product=None, norm=- 2, maxiter=10, rtol=0.001)

computes epsilon-pseudospectrum for the eigenproblem.

Uses the algorithm described in section 5 of

Embree & Keeler (2017). SIAM J. Matrix Anal. Appl. 38, 3: 1028-1054.

to enable the approximation of epsilon-pseudospectra for arbitrary differential-algebraic equation systems.

kint

number of eigenmodes in invariant subspace

zgridtuple

(real, imag) points

mucomplex

center point for pseudospectrum.

pencilint

pencil holding EVP

inner_productfunction

a function that takes two field systems and computes their inner product

compute_mass_matrix(self, Q, inner_product)

Compute the mass matrix M using a given inner product

M must be hermitian, so we compute only half the inner products.

Parameters
  • Q (ndarray) – Matrix of eigenvectors

  • inner_product (function) – a function that takes two field systems and computes their inner product

Returns

Return type

ndarray

set_state(self, system, evector)

Set system to given evector

Parameters
  • system (FieldSystem) – system to fill in

  • evector (ndarray) – eigenvector

plot_spectrum(self, axes=None, spectype='good', xlog=True, ylog=True, real_label='real', imag_label='imag')

Plots the spectrum.

The spectrum plots real parts on the x axis and imaginary parts on the y axis.

Parameters
  • spectype ({‘good’, ‘low’, ‘high’}, optional) – specifies whether to use good, low, or high eigenvalues

  • xlog (bool, optional) – Use symlog on x axis

  • ylog (bool, optional) – Use symlog on y axis

  • real_label (str, optional) – Label to be applied to the real axis

  • imag_label (str, optional) – Label to be applied to the imaginary axis

plot_drift_ratios(self, axes=None)

Plot drift ratios (both ordinal and nearest) vs. mode number.

The drift ratios give a measure of how good a given eigenmode is; this can help set thresholds.

Returns

Return type

matplotlib.figure.Figure