Gallery
NEP-PACK provides a way to access applications in order to easily improve and compare algorithms on realistic problems. Below you will find a description of:
- Standard native gallery. This gallery provides are accessed typically through the function call
nep=nep_gallery(str::String)
and returns a NEP-object with efficient compute-functions. Thestr
is an identifier for the problem. - Berlin-Manchester problem collection. You can access problems from this MATLAB-collection by the call
nep=nep_gallery(NLEVP_NEP,str::String)
, wherestr
is the identified in the MATLAB-package. This requires that you have MATLAB and the problem collection installed. - Extra problems. We provide careful implementations of certain large-scale problems which are too large to include in the main gallery.
julia> nep=nep_gallery("dep0")
julia> λ,v=newton(nep)
(-0.3587189459686265 + 0.0im, Complex{Float64}[0.284742+0.0im, -0.143316+0.0im, 0.278378+0.0im, -0.5009+0.0im, -0.613634+0.0im])
julia> norm(compute_Mlincomb(nep,λ,v))
4.718447854656915e-16
Standard native gallery
nep=nep_gallery(name)
nep=nep_gallery(name,params)
nep=nep_gallery(name,params;kwargs)
Collection of nonlinear eigenvalue problems. Returns a NEP object from a gallery of examples of nonlinear eigenvalue problems. The parameter name
decides which NEP.
Supported problems:
The following list describes the NEP with a certain name
and the associated parameters (params
) and keyword arguments (kwargs
), if any.
dep0
: A randomDEP
with one delay tau = 1, generated with a pseudorandom number generator.
params: size (default = 5)dep0_sparse
: A randomDEP
with sparse matrices and one delay tau = 1, generated with a pseudorandom number generator.
Params: Two optional params determining the size (default = 5) and the fill (default = 0.25)dep0_tridiag
: A randomDEP
with sparse tridiagonal matrices and one delay tau = 1, generated with a pseudorandom number generator.
Params: size (default = 100)dep_symm_double
: ADEP
with double eigenvalues and sparse symmetric matrices and one delay tau = 1.
Params: size (default = 100)
Reference: Example from H. Voss and M. M. Betcke, Restarting iterative projection methods for Hermitian nonlinear eigenvalue problems with minmax property, Numer. Math., 2017dep_double
: ADEP
with a double non-semisimple eigenvalue in λ=3πi.
Reference: Example from E. Jarlebring, Convergence factors of Newton methods for nonlinear eigenvalue problems, LAA, 2012dep1
: ADEP
with one eigenvalue equal to one.pep0
: A randomPEP
, generated with a pseudorandom number generator
Params: size (default = 200)pep0_sym
: A random symmetricPEP
, generated with a pseudorandom number generator
Params: size (default = 200)pep0_sparse
: A randomPEP
with sparse matrices, generated with a pseudorandom number generator
Params: Two optionalparams
determining the size (default = 200) and the fill (default = 0.03)real_quadratic
: A quadraticPEP
with real eigenvalues.
Solutions: Four smallest eigenvalues of the problem:
∘-2051.741417993845
∘-182.101627437811
∘-39.344930222838
∘-4.039879577113
dep_distributed
: A NEP corresponding to a distributed delay time-delay system
Reference: Example in E. Jarlebring and W. Michiels and K. Meerbergen, The infinite Arnoldi method and an application to time-delay systems with distributed delays, Delay Systems - Methods, Applications and New Trends, 2012
Solutions: Some correct eigenvalues:
∘-0.400236388049641 + 0.970633098237807i
∘-0.400236388049641 - 0.970633098237807i
∘2.726146249832675 + 0.000000000000000i
∘-1.955643591177653 + 3.364550574688863i
∘-1.955643591177653 - 3.364550574688863i
∘4.493937056300693 + 0.000000000000000i
∘-1.631513006819252 + 4.555484848248613i
∘-1.631513006819252 - 4.555484848248613i
∘-1.677320660400946 + 7.496870451838560i
∘-1.677320660400946 - 7.496870451838560i
qdep0
: A quadratic delay eigenvalue problem.
Reference: S. W. Gaaf and E. Jarlebring, The infinite Bi-Lanczos method for nonlinear eigenvalue problems, SIAM J. Sci. Comput., 2017qdep1
: A quadratic delay eigenvalue problem.
Reference: E. Jarlebring and W. Michiels and K. Meerbergen, A linear eigenvalue algorithm for the nonlinear eigenvalue problem, Numer. Math., 2011qep_fixed_eig
: A quadratic eigenvalue problem with chosen eigenvalues.
Params: Two optionalparams
determining the size (default = 5) and a vector containing the eigenvalues (default = uniform [-1,1])neuron0
: ADEP
that stems from the modling of a coupled neuron.
Reference: L. P. Shayer and S. A. Campbell, Stability, bifurcation and multistability in a system of two coupled neurons with multiple time delays, SIAM J. Applied Mathematics, 2000. It is also a benchmark example in DDE-BIFTOOL.schrodinger_movebc
: A NEP stemming from the discretization of a Schrödinger equation as described in the NEP-PACK online tutorial. The nonlinearity contains $sinh()$, $cosh()$ and $sqrt()$.
Params: The optional parameters are size of discretizationn
and domain and potential descriptionL0
,L1
,α
andV0
.beam
: ADEP
modelling a beam with delayed stabilizing feedback described. The A1-term is rank one.
Params: size of the matrix (defalut = 100)
Reference: R. Van Beeumen, E. Jarlebring, and W. Michiels, A rank-exploiting infinite Arnoldi algorithm for nonlinear eigenvalue problems, 2016.sine
: A NEP formed by the sum of a polynomial and a sine-function. The sine-term has a rank-one matrix coefficient.
Reference: R. Van Beeumen, E. Jarlebring, and W. Michiels, A rank-exploiting infinite Arnoldi algorithm for nonlinear eigenvalue problems, 2016.bem_fichera
: Represents a boundary element discretization of Helmholtz equation for a domain consisting of the unit cube, except one removed corner (Fichera corner). The mesh is hardcoded. Params: The parameterN
determines the size of the problem (defaultN = 5
).
Reference: The model stems from the model in these papers: A boundary element method for solving PDE eigenvalue problems, M. Steinlechner, bachelor thesis, ETH Zürich, 2010 and Effenberger and Kressner, Chebyshev interpolation for nonlinear eigenvalue problems, BIT Numerical Mathematics, December 2012, Volume 52, Issue 4, pp 933–951dtn_dimer
: NEP with quotients of Bessel functions stemming from the modeling of resonances
Params: This NEP takes two parameters:data_dir::String
andl::Int
. Thedata_dir
specifies the directory of the dowloaded FEM-matrices (available here https://umu.app.box.com/s/b52yux3z9rcl8y0l7la22k0vi062cvu5). The integerl
specifies the number of DtN-terms: 2l+1.
Reference: J. Araujo-Cabarcas, C. Engström and E. Jarlebring, Efficient resonance computations for Helmholtz problems based on a Dirichlet-to-Neumann map, J. Comput. Appl. Math., 330:177-192, 2018)
The MATLAB-package described in T. Betcke, N. J. Higham, V. Mehrmann, Ch. Schröder, F. Tisseur, NLEVP: A Collection of Nonlinear Eigenvalue Problems, ACM Transactions on Mathematical Software 39(2), January 2011 provides a number of benchmark problems for NEPs. These are available in NEP-PACK in two different ways. We have native implementations of some problems (referred to as nlevp_native_
) and the separate GalleryNLEVP
. The native implementation is preferred since the GalleryNLEVP
interfaces with MATLAB and is therefore considerably slower.
nlevp_native_gun
: The benchmark problem from the NLEVP-collection called "gun", represented in the native NEP-PACK format. B.-S. Liao, Z. Bai, L.-Q. Lee, and K. Ko. Nonlinear Rayleigh-Ritz iterative method for solving large scale nonlinear eigenvalue problems. Taiwan. Journal of Mathematics, 14(3):869–883, 2010nlevp_native_cd_player
: The benchmark problem from the NLEVP-collection called "cd_player", represented in the native NEP-PACK format. Y. Chahlaoui, and P. M. Van Dooren, Benchmark examples for model reduction of linear time- invariant dynamical systems. In Dimension Reduction of Large-Scale Systems, P. Benner, V. Mehrmann, and D. C. Sorensen, Eds. Lecture Notes in Computational Science and Engineering Series, vol. 45. Springer-Verlag, Berlin, 380–392, 2005. and P. M. R. Wortelboer, M. Steinbuch, and O. H. Bosgra, Closed-loop balanced reduction with application to a compact disc mechanism. In Selected Topics in Identification, Modeling and Control. Vol. 9. Delft University Press, 47–58, 1996. and W. Draijer, M. Steinbuch, and O. H. Bosgra, Adaptive control of the radial servo system of a compact disc player. Automatica 28, 3, 455–462. 1992.nlevp_native_fiber
: The benchmark problem from the NLEVP-collection called "fiber", represented in the native NEP-PACK format. One of terms in this problem is approximated by interpolation, and may not always coincide with the benchmark. L. Kaufman, Eigenvalue problems in fiber optic design. SIAM J. Matrix Anal. Appl. 28, 1, 105–117, 2006. and X. Huang, Z. Bai, and Y. Su, Nonlinear rank-one modification of the symmetric eigenvalue problem. J. Comput. Math. 28, 2, 218–234, 2010nlevp_native_hadeler
: The benchmark problem from the NLEVP-collection called "hadeler", represented in the native NEP-PACK format. The problem is of the form $M(λ)=(e^λ-1)B+A_0+A_2λ^2$.
Hadeler K. P. 1967. Mehrparametrige und nichtlineare Eigenwertaufgaben. Arch. Rational Mech. Anal. 27, 4, 306–328.nlevp_native_pdde_stability
: The benchmark problem from the NLEVP-collection called "pdde_stability", represented in the native NEP-PACK format. This problem is a quadratic eigenvalue with arbitrary given sizen
. See E. Jarlebring, The Spectrum of Delay-Differential Equations: Numerical Methods, Stability and Perturbation, PhD thesis, TU Braunschweig, Institut Computational Mathematics, Germany, 2008 and H. Fassbender, N. Mackey, D. S. Mackey and C. Schroeder, Structured Polynomial Eigenproblems Related to Time-Delay Systems, ETNA, 2008, vol 31, pp 306-330nlevp_native_loaded_string
: The benchmark problem from the NLEVP-collection called "pdde_stability", represented in the native NEP-PACK format. The parameters are (n,kappa,m) where n is the size, and the NEP is a SPMF with rational terms and the coefficient matrices are rank one modifications of Toeplitz matrices.
S. I. Solov"ev. Preconditioned iterative methods for a class of nonlinear eigenvalue problems. Linear Algebra Appl., 415 (2006), pp.210-229.
Example
julia> nep=nep_gallery("dep0",100);
julia> norm(compute_Mlincomb(nep,1.0+1.0im,ones(size(nep,1))))
57.498446538064954
See also the following galleries:
- GalleryNLEVP
- GalleryWaveguide
Berlin-Manchester collection
If MATLAB and the Berlin-Manchester collection are installed, you can access the Berlin-Manchester collection problems with the GalleryNLEVP. This is a wrapper, which makes MATLAB-calls (via Julias MATLAB interoperability package) for every compute-function call. This can be very inefficient due to overhead. You may want to convert your NEP to a native type, e.g., ChebPEP
.
julia> using GalleryNLEVP
julia> nlevp_path="/home/user/myname/work/src/nlevp";
julia> nep=nep_gallery(NLEVP_NEP,"hadeler",nlevp_path);
julia> λ,v=quasinewton(nep,λ=0.2,logger=1,maxit=20,tol=1e-10);
julia> norm(compute_Mlincomb(nep,λ,v))/norm(v)
9.698206079849311e-11
Problems loaded from the Berlin-Manchester collection are NEP-objects where every call to access a function generates a call to an underlying MATLAB-session. Some problems in the Berlin-Manchester collection have native support in NEP-PACK, i.e., avoiding a MATLAB-access in every call; see nep_gallery
above.
Extra gallery problems
Stand-alone implementation of certain larger problems can be accessed in a similar way to the standard native gallery. A native implementation of a waveguide eigenvalue problem can be accessed as.
julia> using GalleryWaveguide
julia> nep=nep_gallery(WEP,benchmark_problem="TAUSCH");