API
Cim.AbstractContour
Cim.quadpts
Cim.cim
Cim.contr_int_ho
Cim.get_quadpts
Cim.get_quadpts
Cim.get_quadpts
Cim.is_inside
Cim.show_contour!
Cim.show_contour!
Cim.show_eigenvalues!
Cim.show_quadpts!
Cim.AbstractContour
— Typecontours that we plan to support: 1. ellipse 2. circle 3. rectangle (not supported yet)
Cim.quadpts
— Typequadpts
Quadrature points of the trapezoidal rule on the contours.
Properties
N
: the number of the quadrature nodesnodes
: quadrature nodes size of N x 2nodes_prime
: derivative of the parametrization size of N x 2
Cim.cim
— Methodcim(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 contournep
: the nonlinear eigenvalue problemD
: the scale of thenep
l
: the number of columns of random matrixn
: 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.
Cim.contr_int_ho
— Methodcontr_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 matricespbar
: the number of the moments (for p = 0, ..., pbar)
Reference
- Stefan Guttel, Francoise Tisseur, The NEP, Acta Numerica, 2017.
Cim.get_quadpts
— Methodget_quadpts(ctr::circle, num_quadpts::Int64)
Get the quadrature points on the circle ctr
. Here we use trapezoidal rule.
Arguments
ctr
: the contour that we discretizenum_quadpts
: the number of the quadrature nodes
TODO: note that rectangle is different with ellipse and circle
Cim.get_quadpts
— Methodget_quadpts(ctr::ellipse, num_quadpts::Int64)
Get the quadrature points on the ellipse ctr
. Here we use trapezoidal rule.
Arguments
ctr
`: the contour that we discretizenum_quadpts
: the number of the quadrature nodes
TODO: note that rectangle is different with ellipse and circle
Cim.is_inside
— Methodis_inside(λ, ctr::ellipse)
is_inside(λ, ctr::circle)
Check if the eigenvalue λ is in the contour ctr
to avoid spurious eigenvalues.
Cim.show_contour!
— Methodshow_contour!(ax, ctr::circle)
Plot the contour ctr::circle
on Axis ax
by using CairoMakie
.
Cim.show_contour!
— Methodshow_contour!(ax, ctr::ellipse)
Plot the contour ctr::ellipse
on Axis ax
by using CairoMakie
.
Cim.show_eigenvalues!
— Methodshow_eigenvalues(ax, eigvals::AbstractArray)
Plot the eigenvalues on complex plane on Axis ax
by using CairoMakie
.
Cim.show_quadpts!
— Methodshow_quadpts!(ax, pts::quadpts)
Plot the quadrature points pts
on Axis ax
by using CairoMakie
.