EVSL  1.1.0
EigenValues Slicing Library
 All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros
Functions
evsl.h File Reference

This file contains function prototypes and constant definitions for EVSL. More...

#include <complex.h>
#include <stdio.h>
#include "struct.h"
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

int ChebLanNr (double *intv, int maxit, double tol, double *vinit, polparams *pol, int *nevOut, double **lamo, double **Wo, double **reso, FILE *fstats)
 
int ChebLanTr (int lanm, int nev, double *intv, int maxit, double tol, double *vinit, polparams *pol, int *nev2, double **vals, double **W, double **resW, FILE *fstats)
 Chebyshev polynomial filtering Lanczos process [Thick restart version]. More...
 
void set_pol_def (polparams *pol)
 set default values for polparams struct. More...
 
int find_pol (double *intv, polparams *pol)
 Sets the values in pol. More...
 
void free_pol (polparams *pol)
 
int ChebSI (int nev, double *intv, int maxit, double tol, double *vinit, polparams *pol, int *nevo, double **lamo, double **Yo, double **reso, FILE *fstats)
 Chebyshev polynomial filtering Subspace Iteration. More...
 
void SetupPolRec (int n, int max_deg, double tol, double lmin, double lmax, BSolDataPol *data)
 
void SetupPolSqrt (int n, int max_deg, double tol, double lmin, double lmax, BSolDataPol *data)
 
void FreeBSolPolData (BSolDataPol *data)
 
void BSolPol (double *b, double *x, void *data)
 
void extrDiagCsr (csrMat *B, double *d)
 
void diagScalCsr (csrMat *A, double *d)
 
int LanBounds (int msteps, double *v, double *lmin, double *lmax)
 
int LanDos (const int nvec, int msteps, const int npts, double *xdos, double *ydos, double *neig, const double *const intv)
 
int LanDosG (const int nvec, int msteps, const int npts, double *xdos, double *ydos, double *neig, const double *const intv)
 
int LanTrbounds (int lanm, int maxit, double tol, double *vinit, int bndtype, double *lammin, double *lammax, FILE *fstats)
 Lanczos process for eigenvalue bounds [Thick restart version]. More...
 
int scalEigVec (int n, int nev, double *Y, double *sqrtdiag)
 
void set_ratf_def (ratparams *rat)
 Sets default values for ratparams struct. More...
 
int find_ratf (double *intv, ratparams *rat)
 
int set_ratf_solfunc (ratparams *rat, csrMat *A, csrMat *B, SolFuncC *funcs, void **data)
 
void free_rat (ratparams *rat)
 
int RatLanNr (double *intv, int maxit, double tol, double *vinit, ratparams *rat, int *nevOut, double **lamo, double **Wo, double **reso, FILE *fstats)
 Rational filtering Lanczos process [NON-restarted version]. More...
 
int RatLanTr (int lanm, int nev, double *intv, int maxit, double tol, double *vinit, ratparams *rat, int *nev2, double **vals, double **W, double **resW, FILE *fstats)
 RatLanTR filtering Lanczos process [Thick restart version]. More...
 
void csr_copy (csrMat *A, csrMat *B, int allocB)
 copy a csr matrix A into B alloB: 0: will not allocate memory for B (have been alloced outside) 1: will allocate memory for B (same size as A) More...
 
int cooMat_to_csrMat (int cooidx, cooMat *coo, csrMat *csr)
 convert coo to csr More...
 
void free_csr (csrMat *csr)
 memory deallocation for csr matrix More...
 
void free_coo (cooMat *coo)
 memory deallocation for coo matrix More...
 
int matvec (char trans, csrMat *A, double *x, double *y)
 
int SetAMatvec (int n, MVFunc func, void *data)
 Set the user-input matvec routine and the associated data for A. Save them in evsldata. More...
 
int SetBMatvec (int n, MVFunc func, void *data)
 Set the user-input matvec routine and the associated data for B. Save them in evsldata. More...
 
int UnsetAMatvec ()
 
int UnsetBMatvec ()
 
int SetAMatrix (csrMat *A)
 Set the matrix A. More...
 
int SetBMatrix (csrMat *B)
 Set the B matrix. More...
 
int SetBSol (SolFuncR func, void *data)
 Set the solve routine and the associated data for B. More...
 
int SetLTSol (SolFuncR func, void *data)
 Set the solve routine for LT. More...
 
int SetASigmaBSol (ratparams *rat, SolFuncC *func, SolFuncC allf, void **data)
 Set the solve routine and the associated data for A-SIGMA*B if func == NULL, set all functions to be allf. More...
 
int SetStdEig ()
 Set the problem to standard eigenvalue problem. More...
 
int SetGenEig ()
 Set the problem to generalized eigenvalue problem. More...
 
int EVSLStart ()
 Initialize evslData. More...
 
int EVSLFinish ()
 Finish EVSL. More...
 
void SetDiagScal (double *ds)
 Set diagonal scaling matrix D. More...
 
int spslicer (double *sli, double *mu, int Mdeg, double *intv, int n_int, int npts)
 given the dos function defined by mu find a partitioning of sub-interval [a,b] of the spectrum so each subinterval has about the same number of eigenvalues Mdeg = degree.. mu is of length Mdeg+1 [0—> Mdeg] on return [ sli[i],sli[i+1] ] is a subinterval (slice). More...
 
int kpmdos (int Mdeg, int damping, int nvec, double *ab, double *mu, double *ecnt)
 This function computes the coefficients of the density of states in the chebyshev basis. It also returns the estimated number of eigenvalues in the interval given by intv. More...
 
void spslicer2 (double *xi, double *yi, int n_int, int npts, double *sli)
 
double evsl_timer ()
 evsl timer for mac More...
 
void rand_double (int n, double *v)
 
void randn_double (int n, double *v)
 
void sort_double (int n, double *v, int *ind)
 
void linspace (double a, double b, int num, double *arr)
 
void StatsPrint (FILE *fstats)
 
void StatsReset ()
 

Detailed Description

This file contains function prototypes and constant definitions for EVSL.

Definition in file evsl.h.

Function Documentation

void BSolPol ( double *  b,
double *  x,
void *  data 
)

Definition at line 61 of file dos_utils.c.

References _BSolDataPol::cc, _BSolDataPol::dd, _BSolDataPol::deg, _BSolDataPol::mu, pnav(), and _BSolDataPol::wk.

Referenced by main().

Here is the call graph for this function:

Here is the caller graph for this function:

int ChebLanNr ( double *  intv,
int  maxit,
double  tol,
double *  vinit,
polparams pol,
int *  nevOut,
double **  lamo,
double **  Wo,
double **  reso,
FILE *  fstats 
)

Definition at line 60 of file cheblanNr.c.

References _polparams::bar, CGS_DGKS(), CGS_DGKS2(), ChebAv(), DAXPY(), DCOPY(), DDOT(), _polparams::deg, DGEMV(), DNRM2(), DSCAL(), evsl_timer(), evsldata, evslstat, _polparams::gam, _evsldata::ifGenEv, Malloc, min, _evsldata::n, NGS_MAX, orthTol, rand_double(), Realloc, SymmTridEig(), SymmTridEigS(), _evslstat::t_iter, and _evslstat::t_ritz.

Referenced by evsl_cheblannr(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

int ChebLanTr ( int  lanm,
int  nev,
double *  intv,
int  maxit,
double  tol,
double *  vinit,
polparams pol,
int *  nev2,
double **  vals,
double **  W,
double **  resW,
FILE *  fstats 
)

Chebyshev polynomial filtering Lanczos process [Thick restart version].

Parameters
[in]lanmDimension of Krylov subspace [restart dimension]
[in]nevEstimate of number of eigenvalues in the interval – ideally nev == exact number or a little larger. This is not used for testing convergence but it helps make decisions as to when to test convergence ChebLanTr attempts to compute all eigenvalues in the interval and stops only when no more eigenvalyes are left. The convergenve test is a very simple one based on the residual norm for the filtered matrix
[in]intvan array of length 4
[intv[0], intv[1]] is the interval of desired eigenvalues
[intv[2], intv[3]] is the global interval of all eigenvalues
it must contain all eigenvalues of A
[in]maxitmax Num of outer Lanczos iterations (restarts) allowed – Each restart may or use the full lanm lanczos steps or fewer.
[in]toltolerance for convergence. stop when ||res||< tol
[in]vinitinitial vector for Lanczos – [optional]
[in]pola struct containing the parameters of the polynomial. This is set up by a call to find_deg prior to calling chenlanTr
[out]nev2Number of eigenvalues/vectors computed
[out]WA set of eigenvectors [n x nev2 matrix] of unit 2-norm for standard eig prob of unit B-norm for generalized eig prob
[out]valsAssociated eigenvalues [nev2 x 1 vector]
[out]resWAssociated residual norms [nev x 1 vector] 2-norm for standard eig prob B-norm for generalized eig prob
[out]fstatsFile stream which stats are printed to
Returns
Returns 0 on success (or if check_intv() is non-positive), -1 if gamB is outside [-1, 1], and 2 if there are no eigenvalues found.
Warning
memory allocation for W/vals/resW within this function

Definition at line 62 of file cheblanTr.c.

References _polparams::bar, Calloc, CGS_DGKS(), CGS_DGKS2(), ChebAv(), DAXPY(), DCOPY(), DDOT(), DGEMM(), DGEMV(), DNRM2(), DSCAL(), evsl_timer(), evsldata, evslstat, _polparams::gam, _evsldata::ifGenEv, Malloc, max, min, _evsldata::n, NGS_MAX, orthTol, rand_double(), Realloc, SymEigenSolver(), and _evslstat::t_iter.

Referenced by evsl_cheblantr(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

int ChebSI ( int  nev,
double *  intv,
int  maxit,
double  tol,
double *  vinit,
polparams pol,
int *  nevo,
double **  lamo,
double **  Yo,
double **  reso,
FILE *  fstats 
)

Chebyshev polynomial filtering Subspace Iteration.

Parameters
nevEstimate of number of eigenvalues in the interval – ideally nev == exact number or a little larger. ChebSI stops when at least nev eigenvalues are found or when no more candidates are left in interval.
intvAn array of length 4 [intv[0], intv[1]] is the interval of desired eigenvalues [intv[2], intv[3]] is the global interval of all eigenvalues it must contain all eigenvalues of A
maxitMax Num of outer subspace iterations allowed
tolTolerance for convergence. stop when ||res||< tol
vinitNev initial vectors [to be made optional]
polA struct containing the parameters of the polynomial..

Modifies:

Parameters
[out]nevoNumber of eigenvalues/vectors computed
[out]YoA set of eigenvectors [n x nev2 matrix]
[out]lamoAssociated eigenvalues [nev2 x 1 vector]
[out]resoAssociated residual norms [nev x 1 vector]
[out]fstatsFile stream which stats are printed to
Warning
Memory allocation for Yo/lamo/reso within this function

*-----------------— memory for converged Ritz pairs (Y, Lam) */

Definition at line 42 of file chebsi.c.

References _polparams::bar, ChebAv(), DCOPY(), DDOT(), _polparams::deg, DGEMM(), evsl_timer(), evsldata, _polparams::gam, Malloc, _evsldata::n, NBUF, orth(), and SymEigenSolver().

Referenced by main().

Here is the call graph for this function:

Here is the caller graph for this function:

int cooMat_to_csrMat ( int  cooidx,
cooMat coo,
csrMat csr 
)

convert coo to csr

Definition at line 135 of file spmat.c.

References _csrMat::a, csr_resize(), _csrMat::ia, _cooMat::ir, _csrMat::ja, _cooMat::jc, _cooMat::ncols, _cooMat::nnz, _cooMat::nrows, sortrow(), and _cooMat::vv.

Referenced by evsl_coo2csr(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

void csr_copy ( csrMat A,
csrMat B,
int  allocB 
)

copy a csr matrix A into B alloB: 0: will not allocate memory for B (have been alloced outside) 1: will allocate memory for B (same size as A)

Definition at line 106 of file spmat.c.

References _csrMat::a, csr_resize(), _csrMat::ia, _csrMat::ja, _csrMat::ncols, and _csrMat::nrows.

Referenced by main().

Here is the call graph for this function:

Here is the caller graph for this function:

void diagScalCsr ( csrMat A,
double *  d 
)

Definition at line 390 of file spmat.c.

References _csrMat::a, _csrMat::ia, _csrMat::ja, and _csrMat::nrows.

Referenced by main().

Here is the caller graph for this function:

double evsl_timer ( )

evsl timer for mac

evsl timer for mac

Definition at line 14 of file mactime.c.

Referenced by CGS_DGKS(), CGS_DGKS2(), ChebAv(), ChebLanNr(), ChebLanTr(), ChebSI(), EVSLStart(), main(), RatFiltApply(), RatLanNr(), RatLanTr(), SetupASIGMABSolDirect(), SetupBSolDirect(), SymEigenSolver(), SymmTridEig(), SymmTridEigS(), and time_seeder().

Here is the caller graph for this function:

int EVSLFinish ( )

Finish EVSL.

Definition at line 48 of file evsl.c.

References _evsldata::Amv, _evsldata::Bmv, _evsldata::Bsol, and _evsldata::LTsol.

Referenced by evsl_finish(), and main().

Here is the caller graph for this function:

int EVSLStart ( )

Initialize evslData.

Definition at line 27 of file evsl.c.

References _evsldata::Amv, _evsldata::Bmv, _evsldata::Bsol, _evsldata::ds, evsl_timer(), _evsldata::ifGenEv, _evsldata::LTsol, _evsldata::n, and StatsReset().

Referenced by evsl_start(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

void extrDiagCsr ( csrMat B,
double *  d 
)

Definition at line 406 of file spmat.c.

References _csrMat::a, _csrMat::ia, _csrMat::ja, and _csrMat::nrows.

Referenced by main().

Here is the caller graph for this function:

int find_pol ( double *  intv,
polparams pol 
)

Sets the values in pol.

Parameters
intvAn array of length 4, [intv[0, intv[1]] is the interval of desired eigenvalues. [intv[2], intv[3]] is the global interval of all eigenvalues it must contain all eigenvalues of A.
polThe polynomial struct to set the values of.
Warning
Set the following values of pol
mu Expansion coefficients of best polynomial found
gam Site of delta function that is expanded. accurate 'balancing' is done: if p(t) is best approximation to delta function at gam then gam is selected so that p(a)=p(b) - within the tolerance tolBal (set in this function to 1.e-10)
bar If $p(\lambda_i)) \geq $ bar, accept eignevalue as belonging to interval; else reject.
deg degree of polynomial
dd half-width and...
cc ... Center of interval containing all eigenvalues [these are used by ChebAv]

Definition at line 361 of file chebpoly.c.

References _polparams::bar, _polparams::cc, chebxPltd(), chext(), dampcf(), _polparams::damping, _polparams::dd, _polparams::deg, _polparams::gam, _polparams::intvtol, Malloc, max, _polparams::max_deg, min, _polparams::min_deg, _polparams::mu, rootchb(), _polparams::thresh_int, and _polparams::type.

Referenced by evsl_find_pol(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

int find_ratf ( double *  intv,
ratparams rat 
)

Parameters
intv= an array of length 4 [intv[0], intv[1]] is the interval of desired eigenvalues [intv[2], intv[3]] is the global interval of all eigenvalues it must contain all eigenvalues of A

OUT:

Parameters
[out]ratthese are set in rat struct:
omega : expansion coefficients of rational filter
zk : location of the poles used
aa : adjusted left endpoint
bb : adjusted right endpoint
dd : half-width and..
cc : ..center of interval

Definition at line 375 of file ratfilter.c.

References _ratparams::aa, _ratparams::ASIGBsol, _ratparams::bb, _ratparams::beta, contQuad(), Malloc, max, _ratparams::method, min, _ratparams::mulp, _ratparams::num, _ratparams::omega, _ratparams::pow, _ratparams::pw, scaleweigthts(), weights(), and _ratparams::zk.

Referenced by evsl_find_rat(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

void free_coo ( cooMat coo)

memory deallocation for coo matrix

Definition at line 126 of file spmat.c.

References _cooMat::ir, _cooMat::jc, and _cooMat::vv.

Referenced by main().

Here is the caller graph for this function:

void free_csr ( csrMat csr)

memory deallocation for csr matrix

Definition at line 91 of file spmat.c.

References _csrMat::a, _csrMat::ia, _csrMat::ja, and _csrMat::owndata.

Referenced by evsl_free_csr(), main(), and SetupASIGMABSolDirect().

Here is the caller graph for this function:

void free_pol ( polparams pol)

Definition at line 488 of file chebpoly.c.

References _polparams::mu.

Referenced by evsl_free_pol(), and main().

Here is the caller graph for this function:

void free_rat ( ratparams rat)

Definition at line 424 of file ratfilter.c.

References _ratparams::ASIGBsol, _ratparams::mulp, _ratparams::omega, and _ratparams::zk.

Referenced by evsl_free_rat(), and main().

Here is the caller graph for this function:

void FreeBSolPolData ( BSolDataPol data)

Definition at line 53 of file dos_utils.c.

References _BSolDataPol::mu, and _BSolDataPol::wk.

Referenced by main().

Here is the caller graph for this function:

int kpmdos ( int  Mdeg,
int  damping,
int  nvec,
double *  intv,
double *  mu,
double *  ecnt 
)

This function computes the coefficients of the density of states in the chebyshev basis. It also returns the estimated number of eigenvalues in the interval given by intv.


Parameters
Mdegdegree of polynomial to be used.
dampingtype of damping to be used [0=none,1=jackson,2=sigma]
nvecnumber of random vectors to use for sampling
intvan array of length 4
[intv[0] intv[1]] is the interval of desired eigenvalues that must be cut (sliced) into n_int sub-intervals
[intv[2],intv[3]] is the global interval of eigenvalues it must contain all eigenvalues of A
[out]muarray of Chebyshev coefficients
[out]ecntestimated num of eigenvalues in the interval of interest

Definition at line 32 of file spslice.c.

References dampcf(), DDOT(), DNRM2(), DSCAL(), evsldata, _evsldata::ifGenEv, Malloc, max, min, _evsldata::n, PI, rand_double(), and time_seeder().

Referenced by evsl_kpm_spslicer(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

int LanBounds ( int  msteps,
double *  v,
double *  lmin,
double *  lmax 
)

Parameters
mstepsnumber of Lanczos steps
*vinitial vector
[out]*lmin,*lmax[lmin lmax] is the desired interval containing

all eigenvalues of matrix pair (A,B)

Definition at line 21 of file lanbounds.c.

References DAXPY(), DCOPY(), DDOT(), DNRM2(), DSCAL(), evsldata, _evsldata::ifGenEv, Malloc, min, _evsldata::n, orthTol, and SymmTridEig().

Here is the call graph for this function:

int LanDos ( const int  nvec,
int  msteps,
int  npts,
double *  xdos,
double *  ydos,
double *  neig,
const double *const  intv 
)

Computes the density of states (DOS, or spectral density)

Parameters
[in]nvecnumber of sample vectors used
[in]mstepsnumber of Lanczos steps
[in]nptsnumber of sample points used for the DOS curve
[in]*intvStores the two intervals of interest \ intv[0:1] = [a b] = interval where DOS is to be computed intv[2:3] = [lambda_min, lambda_max]\
[out]*xdosLength-npts long vector, x-coordinate points for plotting the DOS. Must be preallocated before calling LanDos
[out]*ydosLength-npts long vector, y-coordinate points for plotting the DOS. Must be preallocated before calling LanDos
[out]neigestimated number of eigenvalues
Note
This only works with the standard eigenvalue problem. Use landosG.c /LanDosG for the generalized eigenvalue problem.

Definition at line 40 of file landos.c.

References Calloc, DAXPY(), DCOPY(), DDOT(), DSCAL(), evsldata, _evsldata::ifGenEv, linspace(), Malloc, max, min, _evsldata::n, orthTol, PI, randn_double(), simpson(), and SymmTridEig().

Referenced by main().

Here is the call graph for this function:

Here is the caller graph for this function:

int LanDosG ( const int  nvec,
const int  msteps,
int  npts,
double *  xdos,
double *  ydos,
double *  neig,
const double *const  intv 
)

Computes the density of states (DOS, or spectral density) using Lanczos algorithm for the generalized eigenvalue problem.

Parameters
[in]nvecnumber of sample vectors used
[in]mstepsnumber of Lanczos steps
[in]nptsnumber of sample points used for the DOS curve
[in]*intvStores the the intervals of interest intv[0:1] = [a b] = interval where DOS is to be computed intv[2:3] = [lambda_min, lambda_max] \
[out]*xdosLength-npts long vector, x-coordinate points for plotting the DOS. Must be preallocated before calling LanDos
[out]*ydosLength-npts long vector, y-coordinate points for plotting the DOS. Must be preallocated before calling LanDos
[out]neigestimated number of eigenvalues
Note
This works for both the standard and generalized eigenvalue problems.

landos.c/LanDos is only for the standard eigenvalue problem.

Definition at line 43 of file landosG.c.

References Calloc, CGS_DGKS(), CGS_DGKS2(), DAXPY(), DCOPY(), DDOT(), DNRM2(), DSCAL(), evsldata, _evsldata::ifGenEv, linspace(), Malloc, max, min, _evsldata::n, NGS_MAX, orthTol, PI, rand_double(), randn_double(), simpson(), and SymmTridEig().

Referenced by main().

Here is the call graph for this function:

Here is the caller graph for this function:

int LanTrbounds ( int  lanm,
int  maxit,
double  tol,
double *  vinit,
int  bndtype,
double *  lammin,
double *  lammax,
FILE *  fstats 
)

Lanczos process for eigenvalue bounds [Thick restart version].

Parameters
[in]lanmDimension of Krylov subspace [restart dimension]
[in]maxitmax Num of outer Lanczos iterations (restarts) allowed – Each restart may or use the full lanm lanczos steps or fewer.
[in]toltolerance for convergence
[in]vinitinitial vector for Lanczos – [optional]
[in]bndtypeType of bound >1 for kato-temple, otherwise simple
[out]lamminLower bound of the spectrum
[out]lammaxUpper bound of the spectrum
[out]fstatsFile stream which stats are printed to
Returns
Returns 0 on success

Definition at line 38 of file lanTrbounds.c.

References Calloc, CGS_DGKS(), CGS_DGKS2(), DAXPY(), DCOPY(), DDOT(), DGEMV(), DNRM2(), DSCAL(), evsldata, _evsldata::ifGenEv, Malloc, min, _evsldata::n, NGS_MAX, orthTol, rand_double(), and SymEigenSolver().

Referenced by evsl_lanbounds(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

void linspace ( double  a,
double  b,
int  num,
double *  arr 
)

Definition at line 54 of file vect.c.

int matvec ( char  trans,
csrMat A,
double *  x,
double *  y 
)
void rand_double ( int  n,
double *  v 
)

Definition at line 11 of file vect.c.

Referenced by ChebLanNr(), ChebLanTr(), evsl_cheblannr(), evsl_cheblantr(), evsl_lanbounds(), evsl_ratlannr(), evsl_ratlantr(), kpmdos(), LanDosG(), LanTrbounds(), main(), RatLanNr(), and RatLanTr().

Here is the caller graph for this function:

void randn_double ( int  n,
double *  v 
)

Definition at line 26 of file vect.c.

Referenced by LanDos(), and LanDosG().

Here is the caller graph for this function:

int RatLanNr ( double *  intv,
int  maxit,
double  tol,
double *  vinit,
ratparams rat,
int *  nevOut,
double **  lamo,
double **  Wo,
double **  reso,
FILE *  fstats 
)

Rational filtering Lanczos process [NON-restarted version].


Parameters
intvan array of length 4 [intv[0], intv[1]] is the interval of desired eigenvalues [intv[2], intv[3]] is the global interval of all eigenvalues it must contain all eigenvalues of A
maxitmax Num of outer Lanczos steps allowed –[max dim of Krylov subspace]
toltolerance for convergence. The code uses a stopping criterion based on the convergence of the restricted trace. i.e., the sum of the eigenvalues of T_k that are in the desired interval. This test is rather simple since these eigenvalues are above `bar'. We want the relative error on this restricted trace to be less than tol. Note that the test performed on filtered matrix only - but the actual residual norm associated with the original matrix A is returned
vinitinitial vector for Lanczos – [optional]
Warning
RatLanNr() Modifies the following variables:
Parameters
[out]ratA struct containing the rational filter
[out]nevOutNumber of eigenvalues/vectors computed
[out]WoA set of eigenvectors [n x nevOut matrix]
[out]lamoAssociated eigenvalues [nevOut x 1 vector]
[out]resoAssociated residual norms [nev x 1 vector]
[out]fstatsFile stream which stats are printed to

Definition at line 54 of file ratlanNr.c.

References CGS_DGKS(), CGS_DGKS2(), DAXPY(), DCOPY(), DDOT(), DGEMV(), DNRM2(), _evsldata::ds, DSCAL(), evsl_timer(), evsldata, evslstat, _evsldata::ifGenEv, Malloc, min, _evsldata::n, NGS_MAX, orthTol, rand_double(), RatFiltApply(), Realloc, SymmTridEig(), SymmTridEigS(), _evslstat::t_iter, and _evslstat::t_ritz.

Referenced by evsl_ratlannr(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

int RatLanTr ( int  lanm,
int  nev,
double *  intv,
int  maxit,
double  tol,
double *  vinit,
ratparams rat,
int *  nev2,
double **  vals,
double **  W,
double **  resW,
FILE *  fstats 
)

RatLanTR filtering Lanczos process [Thick restart version].

Parameters
[in]lanmDimension of Krylov subspace [restart dimension]
[in]nevEstimate of number of eigenvalues in the interval – ideally nev == exact number or a little larger. This is not used for testing convergence but it helps make decisions as to when to test convergence RatLanTr attempts to compute all eigenvalues in the interval and stops only when no more eigenvalyes are left. The convergenve test is a very simple one based on the residual norm for the filtered matrix
[in]intvan array of length 4
[intv[0], intv[1]] is the interval of desired eigenvalues
[intv[2], intv[3]] is the global interval of all eigenvalues
it must contain all eigenvalues of A
[in]maxitmax Num of outer Lanczos iterations (restarts) allowed – Each restart may or use the full lanm lanczos steps or fewer.
[in]toltolerance for convergence. stop when ||res||< tol
[in]vinitinitial vector for Lanczos – [optional]
[in]rata struct containing the parameters of the rational.
[out]nev2Number of eigenvalues/vectors computed
[out]WA set of eigenvectors [n x nev2 matrix] of unit 2-norm for standard eig prob of unit B-norm for generalized eig prob
[out]valsAssociated eigenvalues [nev2 x 1 vector]
[out]resWAssociated residual norms [nev x 1 vector] 2-norm for standard eig prob B-norm for generalized eig prob
[out]fstatsFile stream which stats are printed to
Returns
Returns 0 on success (or if check_intv() is non-positive), and 2 if there are no eigenvalues found.
Warning
memory allocation for W/vals/resW within this function

Definition at line 61 of file ratlanTr.c.

References Calloc, CGS_DGKS(), CGS_DGKS2(), DAXPY(), DCOPY(), DDOT(), DGEMM(), DGEMV(), DNRM2(), _evsldata::ds, DSCAL(), evsl_timer(), evsldata, evslstat, _evsldata::ifGenEv, Malloc, max, min, _evsldata::n, NGS_MAX, orthTol, rand_double(), RatFiltApply(), Realloc, SymEigenSolver(), and _evslstat::t_iter.

Referenced by evsl_ratlantr(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

int scalEigVec ( int  n,
int  nev,
double *  Y,
double *  sqrtdiag 
)
void set_pol_def ( polparams pol)

set default values for polparams struct.

Definition at line 18 of file chebpoly.c.

References _polparams::damping, _polparams::deg, _polparams::intvtol, _polparams::max_deg, _polparams::min_deg, _polparams::mu, _polparams::thresh_ext, _polparams::thresh_int, and _polparams::tol.

Referenced by evsl_find_pol(), and main().

Here is the caller graph for this function:

void set_ratf_def ( ratparams rat)

Sets default values for ratparams struct.

Definition at line 345 of file ratfilter.c.

References _ratparams::aa, _ratparams::bar, _ratparams::bb, _ratparams::beta, _ratparams::method, _ratparams::num, and _ratparams::pw.

Referenced by evsl_find_rat(), and main().

Here is the caller graph for this function:

int set_ratf_solfunc ( ratparams rat,
csrMat A,
csrMat B,
SolFuncC funcs,
void **  data 
)
int SetAMatrix ( csrMat A)

Set the matrix A.

Definition at line 68 of file evsl.c.

References _evsldata::Amv, Calloc, _EVSLMatvec::data, _EVSLMatvec::func, matvec_csr(), _evsldata::n, and _csrMat::ncols.

Referenced by evsl_seta_csr(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

int SetAMatvec ( int  n,
MVFunc  func,
void *  data 
)

Set the user-input matvec routine and the associated data for A. Save them in evsldata.

Warning
Once this matvec func is set, matrix A will be ignored even it is provided

Definition at line 100 of file evsl.c.

References _evsldata::Amv, Calloc, _EVSLMatvec::data, _EVSLMatvec::func, and _evsldata::n.

Referenced by evsl_setamv(), and main().

Here is the caller graph for this function:

int SetASigmaBSol ( ratparams rat,
SolFuncC func,
SolFuncC  allf,
void **  data 
)

Set the solve routine and the associated data for A-SIGMA*B if func == NULL, set all functions to be allf.

Definition at line 168 of file evsl.c.

References _ratparams::ASIGBsol, Calloc, _EVSLASIGMABSol::data, _EVSLASIGMABSol::func, and _ratparams::num.

Referenced by main(), and set_asigmabsol_direct().

Here is the caller graph for this function:

int SetBMatrix ( csrMat B)

Set the B matrix.

Definition at line 83 of file evsl.c.

References _evsldata::Bmv, Calloc, _EVSLMatvec::data, _EVSLMatvec::func, matvec_csr(), _evsldata::n, and _csrMat::ncols.

Referenced by evsl_setb_csr(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

int SetBMatvec ( int  n,
MVFunc  func,
void *  data 
)

Set the user-input matvec routine and the associated data for B. Save them in evsldata.

Warning
Once this matvec func is set, matrix B will be ignored even it is provided

Definition at line 117 of file evsl.c.

References _evsldata::Bmv, Calloc, _EVSLMatvec::data, _EVSLMatvec::func, and _evsldata::n.

Referenced by evsl_setbmv().

Here is the caller graph for this function:

int SetBSol ( SolFuncR  func,
void *  data 
)

Set the solve routine and the associated data for B.

Definition at line 132 of file evsl.c.

References _evsldata::Bsol, Calloc, _EVSLBSol::data, and _EVSLBSol::func.

Referenced by evsl_setbsol(), main(), and setup_bsol_direct().

Here is the caller graph for this function:

void SetDiagScal ( double *  ds)

Set diagonal scaling matrix D.

Definition at line 197 of file evsl.c.

References _evsldata::ds.

int SetGenEig ( )

Set the problem to generalized eigenvalue problem.

Definition at line 158 of file evsl.c.

References _evsldata::ifGenEv.

Referenced by evsl_set_geneig(), and main().

Here is the caller graph for this function:

int SetLTSol ( SolFuncR  func,
void *  data 
)

Set the solve routine for LT.

Definition at line 185 of file evsl.c.

References Calloc, _EVSLLTSol::data, _EVSLLTSol::func, and _evsldata::LTsol.

Referenced by evsl_setltsol(), main(), and setup_bsol_direct().

Here is the caller graph for this function:

int SetStdEig ( )

Set the problem to standard eigenvalue problem.

Definition at line 148 of file evsl.c.

References _evsldata::ifGenEv.

Referenced by main().

Here is the caller graph for this function:

void SetupPolRec ( int  n,
int  max_deg,
double  tol,
double  lmin,
double  lmax,
BSolDataPol data 
)

Definition at line 37 of file dos_utils.c.

References rec(), and SetupBPol().

Referenced by main().

Here is the call graph for this function:

Here is the caller graph for this function:

void SetupPolSqrt ( int  n,
int  max_deg,
double  tol,
double  lmin,
double  lmax,
BSolDataPol data 
)

Definition at line 45 of file dos_utils.c.

References isqrt(), and SetupBPol().

Referenced by main().

Here is the call graph for this function:

Here is the caller graph for this function:

void sort_double ( int  n,
double *  v,
int *  ind 
)

Definition at line 92 of file vect.c.

References compare1(), compare2(), _doubleint::d, _doubleint::i, and Malloc.

Referenced by exeiglap3(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

int spslicer ( double *  sli,
double *  mu,
int  Mdeg,
double *  intv,
int  n_int,
int  npts 
)

given the dos function defined by mu find a partitioning of sub-interval [a,b] of the spectrum so each subinterval has about the same number of eigenvalues Mdeg = degree.. mu is of length Mdeg+1 [0—> Mdeg] on return [ sli[i],sli[i+1] ] is a subinterval (slice).


Parameters
*slisee above (output)
*mucoeffs of polynomial (input)
Mdegdegree of polynomial (input)
*intvan array of length 4 [intv[0] intv[1]] is the interval of desired eigenvalues that must be cut (sliced) into n_int sub-intervals [intv[2],intv[3]] is the global interval of eigenvalues it must contain all eigenvalues of A
n_intnumber of slices wanted (input)
nptsnumber of points to use for discretizing the interval [a b]. The more points the more accurate the intervals. it is recommended to set npts to a few times the number of eigenvalues in the interval [a b] (input).

Definition at line 204 of file spslice.c.

References intChx(), linspace(), Malloc, max, min, and save_vec().

Referenced by evsl_kpm_spslicer(), and main().

Here is the call graph for this function:

Here is the caller graph for this function:

void spslicer2 ( double *  xi,
double *  yi,
int  n_int,
int  npts,
double *  sli 
)

Interval partitioner based for Lanczos DOS output

Parameters
[in]xicoordinates of interval [a b]
[in]yiyi[k] = integral of the does from a to xi[k]
[in]n_intNumber of desired sub-intervals
[in]nptsnumber of integration points (length of xi)
[out]sliArray of length n_int containing the boundaries of the intervals. [sli[i], [sli[i+1]] is the i-th interval with sli[0] = xi[0] and sli[n_int] = xi[npts-1]

Definition at line 29 of file spslice2.c.

References simpson().

Referenced by main().

Here is the call graph for this function:

Here is the caller graph for this function:

void StatsPrint ( FILE *  fstats)
void StatsReset ( )

Definition at line 79 of file stats.c.

References evslstat.

Referenced by EVSLStart(), and main().

Here is the caller graph for this function:

int UnsetAMatvec ( )
int UnsetBMatvec ( )