API

Contour

Quadrature

Cim.quadptsType
quadpts

Quadrature points of the trapezoidal rule on the contours.

Properties

  • N: the number of the quadrature nodes
  • nodes: quadrature nodes size of N x 2
  • nodes_prime: derivative of the parametrization size of N x 2
source
Cim.get_quadptsFunction
get_quadpts(ctr::ellipse, num_quadpts::Int64) 
get_quadpts(ctr::circle, num_quadpts::Int64)

Get the quadrature points on the contour (ellipse or circle) ctr. Here we use the composite trapezoidal rule.

Arguments

  • ctr: the contour that we discretize
  • num_quadpts: the number of the quadrature nodes
source

Nonlinear eigenvalue problems

Cim.QepType
Qep{T}

Construct quadratic eigenvlaue problems

$(\mathbf{A}_{0} + \alpha \mathbf{A}_{1} + \alpha^2 \mathbf{A}_{2}) \mathbf{u} = 0.$

Fields

  • A₀: the matrix in the zero order term
  • A₁: the matrix in the first order term
  • A₂: the matrix in the second order term
source

Contour integral method

Cim.cimFunction
cim(ctr::AbstractContour, nep::Function, D, l::Int64; n=50, tol=1e-12)

Contour integral method to calculate the eigenvalues inside the contour ctr

Arugments

  • ctr: the contour
  • nep: the nonlinear eigenvalue problem
  • d: the scale of the nep
  • l: the number of columns of random matrix
  • n: the number of the quadrature points (here we use the trapezoid rule)
  • tol: tolerance to determine the number of nonzero singular values

Reference

  • Wolf-Jurgen Beyn, An integral method for solving NEPs, Linear Algebra Appl., 2012.
source
cim(ctr::AbstractContour, nep::Qep, d::Int64, l::Int64; n=50, tol=1e-12)

Use the contour integral method to solve quadratic eigenvalue problems.

source
Cim.contr_int_hoFunction
contr_int_ho()

Contour integral method with high-order moments

Arguments

In this function, we follow the notations in Stefan Guttel, Francoise Tisseur, Acta Numerica, 2017

  • ctr:
  • nep:
  • D:
  • r, l: size of the probing matrices
  • pbar: the number of the moments (for p = 0, ..., pbar)

Reference

  • Stefan Guttel, Francoise Tisseur, The NEP, Acta Numerica, 2017.
source

Visualization

Cim.show_contour!Function
show_contour!(ax, ctr::ellipse)

Plot the contour ctr::ellipse on Axis ax by using CairoMakie.

source
show_contour!(ax, ctr::circle)

Plot the contour ctr::circle on Axis ax by using CairoMakie.

source
Cim.show_eigenvalues!Function
show_eigenvalues(ax, eigvals::AbstractArray)

Plot the eigenvalues on complex plane on Axis ax by using CairoMakie.

source
Cim.show_quadpts!Function
show_quadpts!(ax, pts::quadpts)

Plot the quadrature points pts on Axis ax by using CairoMakie.

source

Utils

Cim.is_insideFunction
is_inside(λ, ctr::ellipse)
is_inside(λ, ctr::circle)

Check if the eigenvalue λ is in the contour ctr to avoid spurious eigenvalues.

source
Cim.print_vecFunction
print_vec(vec; precision::Int=6)

Print the vector vec to the REPL.

source