naginterfaces.library.sparseig.feast_​complex_​herm_​solve

naginterfaces.library.sparseig.feast_complex_herm_solve(handle, irevcm, ze, x, y, m0, d, z, eps, itera, resid, io_manager=None)[source]

feast_complex_herm_solve is an iterative solver used to find some of the eigenvalues and the corresponding eigenvectors of a standard or generalized eigenvalue problem defined by complex Hermitian matrices. This is part of a suite of functions that also includes feast_init(), feast_option() and feast_symm_contour().

For full information please refer to the NAG Library document for f12jr

https://support.nag.com/numeric/nl/nagdoc_30.2/flhtml/f12/f12jrf.html

Parameters
handleHandle

The handle to the internal data structure used by the NAG FEAST suite. It needs to be initialized by feast_init(). It must not be changed between calls to the NAG FEAST suite.

irevcmint

On initial entry: , otherwise an error condition will be raised.

On intermediate entry: must be unchanged from its previous exit value. Changing to any other value between calls will result in an error.

Note: the matrices , and referred to in this section are all of size and are stored in the arrays , and , respectively.

zecomplex

On initial entry: need not be set.

xcomplex, ndarray, shape , modified in place

On initial entry: need not be set.

On intermediate exit:

if , the calling program must compute , storing the result in prior to re-entry.

If , the calling program must compute , storing the result in prior to re-entry.

Note: the matrices and are stored in the first columns of the arrays and , respectively.

ycomplex, ndarray, shape , modified in place

On initial entry: need not be set.

On intermediate exit:

if , the calling program must compute the solution to the linear system , overwriting with the result , prior to re-entry. The linear system has right-hand sides.

If , the calling program must compute the solution to the linear system , overwriting with the result , prior to re-entry. The linear system has right-hand sides.

m0int

On initial entry: the size of the search subspace used to find the eigenvalues. This should exceed the number of eigenvalues within the search contour. See Further Comments for further details.

On intermediate entry: must remain unchanged.

dfloat, ndarray, shape , modified in place

On initial entry: if the option was set using the option setting function feast_option(), then should contain an initial guess at the eigenvalues lying within the eigenvector search subspace (this subspace should be specified by ), otherwise need not be set.

On final exit: the first entries in contain the eigenvalues.

Note: if the option was set using the option setting function feast_option(), then on final exit contains an estimate of the eigenvalues after a single contour integral.

zcomplex, ndarray, shape , modified in place

On initial entry: if the option was set using the option setting function feast_option(), then should contain an initial guess at the eigenvector search subspace, otherwise need not be set.

On intermediate exit: must not be changed.

On final exit: the first columns of contain the eigenvectors corresponding to the eigenvalues found within the contour.

Note: if the option was set using the option setting function feast_option(), then on final exit columns of contain the current search subspace after one contour integral.

epsfloat

On initial entry: need not be set.

iteraint

On initial entry: need not be set.

residfloat, ndarray, shape , modified in place

On initial entry: need not be set.

On final exit: for , contains the relative residual, in the -norm, of the th eigenpair found, that is , where and are the lower and upper bounds respectively of the eigenvalue search interval.

io_managerFileObjManager, optional

Manager for I/O in this routine.

Returns
irevcmint

On intermediate exit: has the following meanings.

The calling program must compute a factorization of the matrix suitable for solving a linear system, for example using sparse.complex_gen_precon_ilu which computes an incomplete factorization of a complex sparse matrix. All arguments to the routine must remain unchanged.

Note: the factorization can be computed in single precision.

The calling program must compute the solution to the linear system , overwriting with the result . The matrix has previously been factorized (when was returned) and this factorization can be reused here.

Note: the solve can be performed in single precision.

Optionally, the calling program must compute a factorization of the matrix . This need only be done if it is not possible to use the factorization computed when was returned to solve linear systems involving .

Note: the solve can be performed in single precision.

The calling program must compute the solution to the linear system , overwriting with the result . If it is not possible to use the factorization of (computed when was returned) then the factorization of (computed when was returned) should be used here.

The calling program must compute , storing the result in .

The calling program must compute , storing the result in . If a standard eigenproblem is being solved (so that ) then the calling program should set .

On final exit: : feast_complex_herm_solve has completed its tasks. The value of determines whether the iteration has been successfully completed, or whether errors have been detected.

Note: the matrices , and referred to in this section are all of size and are stored in the arrays , and , respectively.

zecomplex

On intermediate exit: contains the current point on the contour.

If , then this must be used by the calling program to form a factorization of the matrix .

If , then, optionally, this can be used to form a factorization of .

m0int

If the initial search subspace was found by feast_complex_herm_solve to be too large, then a new smaller suitable choice is returned.

nconvint

The number of eigenvalues found within the search contour.

Note: if the option was set in the option setting function feast_option(), then contains a stochastic estimate of the number of eigenvalues within the search contour.

epsfloat

The relative error on the trace. At iteration , is given by the expression , where is the sum of the eigenvalues found at the th iteration, and and are the lower and upper bounds respectively of the eigenvalue search interval.

iteraint

The number of subspace iterations performed.

Raises
NagValueError
(errno )

Either the contour setting function feast_symm_contour() has not been called prior to the first call of this function or the supplied has become corrupted.

(errno )

An internal error occurred in the reduced eigenvalue solver. A possible cause is that the matrix is not positive definite.

(errno )

On entry, .

Constraint: .

(errno )

On entry, and .

Constraint: .

(errno )

On initial entry, .

Constraint: .

(errno )

On intermediate entry, .

Constraint: , , , , or .

(errno )

The option was set using the option setting function feast_option() but no nonzero elements were found in the supplied subspace.

Warns
NagAlgorithmicWarning
(errno )

No eigenvalues were found within the search contour.

(errno )

The function did not converge after the maximum number of iterations.

(errno )

The size of the eigenvector search subspace, , is too small.

(errno )

The option was set using feast_option(). Columns of contain the search subspace after one contour integral and contains an estimate of the eigenvalues.

(errno )

The option was set using feast_option() so only a stochastic estimate of the number of eigenvalues within the contour has been returned.

Notes

The suite of functions is designed to calculate some of the eigenvalues, , and the corresponding eigenvectors, , of a standard eigenvalue problem , or of a generalized eigenvalue problem , where the coefficient matrices and are sparse Hermitian, and is positive definite. The suite can also be used to find selected eigenvalues/eigenvectors of smaller scale, dense, Hermitian problems.

feast_complex_herm_solve is a reverse communication function, based on the FEAST eigensolver, described in Polizzi (2009), which finds eigenvalues using contour integration. Prior to calling feast_complex_herm_solve, the contour definition function feast_symm_contour() is used to specify a search interval on the real line, , within which eigenvalues will be sought (note that the eigenvalues of complex Hermitian eigenproblems are themselves real). feast_symm_contour() uses this interval to define nodes and weights for an elliptical contour to be used by feast_complex_herm_solve.

The setup function feast_init() and the contour definition function feast_symm_contour() must be called before feast_complex_herm_solve. Between the calls to feast_init() and feast_symm_contour(), options may be set by calls to the option setting function feast_option().

feast_complex_herm_solve uses reverse communication, i.e., it returns repeatedly to the calling program with the argument (see Parameters) set to specified values which require the calling program to carry out one of the following tasks:

  • compute a factorization of the matrix , where is a point on the search contour;

  • optionally, compute a factorization of the matrix (this need only be done if the factorization does not allow linear systems involving to be solved);

  • solve a linear system involving or , using the factorizations above;

  • compute the matrix product ;

  • compute the matrix product ;

  • notify the completion of the computation.

The number of contour points, the number of iterations, and other options can all be set using the option setting function feast_option() (see Other Parameters for feast_option for details on setting options and of the default settings). The search contour itself is defined by a call to feast_symm_contour().

References

Polizzi, E, 2009, Density-Matrix-Based Algorithms for Solving Eigenvalue Problems, Phys. Rev. B. (79), 115112