Deflation

Deflation

Several NEP-algorithms are able to find one eigenvalue, but may have difficulties finding several eigenvalues. Deflation is a transformation technique which can transform a NEP by effectively removing computed eigenvalues, and allowing several eigenvalues to be computed by repeated application of the same NEP-algorithm.

NEP-PACK provides a solver-independent implementation of deflation which can be combined (essentially) with any NEP-solver. NEP-PACK also has some NEP-solver deflation techniques and reconvergence avoidance techniques incoprorated directly in the solver, e.g., in the nonlinear Arnoldi method (nlar), the Jacobi-Davidson method (jd_betcke) and Broyden's method (broyden).

The technique takes a NEP and a solution and creates a bigger NEP with one dimension larger but where the eigenvalue is removed from the solution set. Due to the abstraction of NEP-objects in NEP-PACK, the deflated NEP is again a NEP and we can apply the NEP-solver to the deflated NEP.

Note

More elaborate deflation examples can be found in the tutorial on deflation.

Theory

The theory follows the presentation of the technique in the PhD thesis of Cedric Effenberger. It can be summarized as follows, in a somewhat simplified form (for the index one case). The deflation is based on a theory for NEPs essentially stating that if $(s,x)$ is an eigenpair, then under certain general conditions (which we implicitly assume are satisfied), the extended nonlinear eigenvalue problem

\[T(λ):=\begin{bmatrix}M(λ)&M(λ)x(s-λ)^{-1}\\ x^T & 0\end{bmatrix}\]

has the same eigenvalues as the original problem except for the eigenvalue $s$ which is no longer part of the solution set. We have effectively removed (i.e. deflated) the eigenpair (s,x). More eigenpairs can be deflated with techniques of partial Schur factorizations, which the user does not need to be aware of, due to the abstraction provided by the functions below. When we create a deflated NEP, we create the NEP $T$.

There are several ways to represent the $T$, which is why deflation has several modes. If you run

julia> dnep=deflate_eigpair(nep,λ1,v1,mode=:SPMF)

the dnep will be of the type AbstractSPMF. More precisely, if

\[M(λ)=A_1f_1(λ)+\cdots+A_mf_m(λ)\]

the deflated NEP will be

\[T(λ)= \begin{bmatrix}A_1&0\\0 & 0\end{bmatrix}f_1(λ)+\cdots+ \begin{bmatrix}A_m&0\\0 & 0\end{bmatrix}f_m(λ)+\]
\[\begin{bmatrix}0&A_1x\\0 & 0\end{bmatrix}\frac{f_1(λ)}{s-λ}+\cdots+ \begin{bmatrix}0&A_mx\\0 & 0\end{bmatrix}\frac{f_m(λ)}{s-λ}+ \begin{bmatrix}0&0\\x^T & 0\end{bmatrix}\]

Clearly, the deflated NEP has more SPMF-terms than the original NEP. When the parameter mode=:SPMF is set, the deflation method will explicitly construct an SPMF_NEP. This is not recommended if you have many SPMF-terms in the original problem, but can be efficient when you only have a few terms. (Some additional exploitation is however implemented, since we can use the fact that the introduced terms are of low rank, and therefore naturally represented as a LowRankFactorizedNEP.)

If you select mode=:Generic, the compute functions are implemented without the use of SPMF, and can be more efficient if the NEP has many SPMF-terms. When mode=:MM the compute-functions are all implemented by calls to compute_MM. This will not be efficient if compute_Mder(nep,λ,der) where der>0 is needed.

Functions

dnep=deflate_eigpair(orgnep::NEP,λ,v,[mode=:Auto])

This function creates a deflated NEP based on $(λ,v)$, which are assumed to an eigenpair of nep. Effectively, the function will return dnep::NEP which has the same solutions as orgnep, except those corresponding to $(λ,v)$. Deflation is typically used to avoid reconvergence.

If orgnep is a DeflatedNEP, the orgnep the deflation in orgnep will be updated.

The mode kwarg can be :Auto, :Generic, :SPMF, :MM. This specifies how the deflated NEP should be represented. Which mode is the most efficient depends on many problem properties. If the original NEP is an AbstractSPMF with only a few terms, mode=:SPMF may be efficient. The SPMF-mode is based on a diagonalization of the deflated invariant pair and is not necessarily robust when you deflate eigenvalues near to each other. When mode=:MM is used, all compute functions are implemented via calls to the compute_MM. This can work well for small dense problems. The :Generic is based on an explicit derivation of the problem (via binomial expansion) which can be efficient if low order derivates are needed. If :Auto is selected, NEP-PACK tries to determine which one is the most efficient based on the orgnep.

Example:

julia> nep=nep_gallery("dep0");
julia> (λ,v)=newton(nep,v=ones(size(nep,1)));
julia> dnep=deflate_eigpair(nep,λ,v)
julia> (λ2,v2)=augnewton(dnep,v=ones(size(dnep,1)));  # this converges to different eigval
julia> using LinearAlgebra;
julia> minimum(svdvals(compute_Mder(nep,λ2)))
2.5161012836775824e-17

The function get_deflated_eigpairs() extracts the eigenpairs that have been deflated. The returned pairs are eigenpairs of the original NEP:

julia> dnep=deflate_eigpair(dnep,λ2,v2);
julia> (D,V)=get_deflated_eigpairs(dnep)
julia> norm(compute_Mlincomb(nep,D[1],V[:,1]))
2.3970746442479104e-16
julia> norm(compute_Mlincomb(nep,D[2],V[:,2]))
8.101585801848734e-16

References

  • C. Effenberger, Robust solution methods for nonlinear eigenvalue problems, PhD thesis, 2013, EPF Lausanne
(D,V)=get_deflated_eigpairs(S,V,n)
(D,V)=get_deflated_eigpairs(dnep::DeflatedNEP [λ,v])

Returns a vector of eigenvalues D and a matrix with corresponding eigenvectors V of the invariant pair S,V. The eigenpairs correspond to the original problem, underlying the DeflatedNEP. The optional parameters λ,v allows the inclusion of an additional eigpair. Essentially, the optional parameters are the expanding the deflation and the running get_deflated_eigpairs with kwarg, i.e.,

(D,V)=get_deflated_eigpairs(deflate_eigpair(dnep,λ,v))`

See example in deflate_eigpair.

To the top