EVSL
1.1.0
EigenValues Slicing Library
|
This file contains function prototypes and constant definitions for EVSL. More...
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 () |
This file contains function prototypes and constant definitions for EVSL.
Definition in file evsl.h.
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().
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().
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].
[in] | lanm | Dimension of Krylov subspace [restart dimension] |
[in] | nev | Estimate 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] | 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 |
[in] | maxit | max Num of outer Lanczos iterations (restarts) allowed – Each restart may or use the full lanm lanczos steps or fewer. |
[in] | tol | tolerance for convergence. stop when ||res||< tol |
[in] | vinit | initial vector for Lanczos – [optional] |
[in] | pol | a struct containing the parameters of the polynomial. This is set up by a call to find_deg prior to calling chenlanTr |
[out] | nev2 | Number of eigenvalues/vectors computed |
[out] | W | A set of eigenvectors [n x nev2 matrix] of unit 2-norm for standard eig prob of unit B-norm for generalized eig prob |
[out] | vals | Associated eigenvalues [nev2 x 1 vector] |
[out] | resW | Associated residual norms [nev x 1 vector] 2-norm for standard eig prob B-norm for generalized eig prob |
[out] | fstats | File stream which stats are printed to |
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().
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.
nev | Estimate 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. |
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 |
maxit | Max Num of outer subspace iterations allowed |
tol | Tolerance for convergence. stop when ||res||< tol |
vinit | Nev initial vectors [to be made optional] |
pol | A struct containing the parameters of the polynomial.. |
Modifies:
[out] | nevo | Number of eigenvalues/vectors computed |
[out] | Yo | A set of eigenvectors [n x nev2 matrix] |
[out] | lamo | Associated eigenvalues [nev2 x 1 vector] |
[out] | reso | Associated residual norms [nev x 1 vector] |
[out] | fstats | File stream which stats are printed to |
*-----------------— 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().
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().
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().
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().
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().
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().
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().
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().
int find_pol | ( | double * | intv, |
polparams * | pol | ||
) |
Sets the values in pol.
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. |
pol | The polynomial struct to set the values of. |
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().
int find_ratf | ( | double * | intv, |
ratparams * | rat | ||
) |
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:
[out] | rat | these 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().
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().
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().
void free_pol | ( | polparams * | pol | ) |
Definition at line 488 of file chebpoly.c.
References _polparams::mu.
Referenced by evsl_free_pol(), and main().
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().
void FreeBSolPolData | ( | BSolDataPol * | data | ) |
Definition at line 53 of file dos_utils.c.
References _BSolDataPol::mu, and _BSolDataPol::wk.
Referenced by main().
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.
Mdeg | degree of polynomial to be used. | |
damping | type of damping to be used [0=none,1=jackson,2=sigma] | |
nvec | number of random vectors to use for sampling | |
intv | an 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] | mu | array of Chebyshev coefficients |
[out] | ecnt | estimated 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().
int LanBounds | ( | int | msteps, |
double * | v, | ||
double * | lmin, | ||
double * | lmax | ||
) |
msteps | number of Lanczos steps | |
*v | initial 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().
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)
[in] | nvec | number of sample vectors used |
[in] | msteps | number of Lanczos steps |
[in] | npts | number of sample points used for the DOS curve |
[in] | *intv | Stores 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] | *xdos | Length-npts long vector, x-coordinate points for plotting the DOS. Must be preallocated before calling LanDos |
[out] | *ydos | Length-npts long vector, y-coordinate points for plotting the DOS. Must be preallocated before calling LanDos |
[out] | neig | estimated number of eigenvalues |
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().
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.
[in] | nvec | number of sample vectors used |
[in] | msteps | number of Lanczos steps |
[in] | npts | number of sample points used for the DOS curve |
[in] | *intv | Stores 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] | *xdos | Length-npts long vector, x-coordinate points for plotting the DOS. Must be preallocated before calling LanDos |
[out] | *ydos | Length-npts long vector, y-coordinate points for plotting the DOS. Must be preallocated before calling LanDos |
[out] | neig | estimated number of eigenvalues |
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().
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].
[in] | lanm | Dimension of Krylov subspace [restart dimension] |
[in] | maxit | max Num of outer Lanczos iterations (restarts) allowed – Each restart may or use the full lanm lanczos steps or fewer. |
[in] | tol | tolerance for convergence |
[in] | vinit | initial vector for Lanczos – [optional] |
[in] | bndtype | Type of bound >1 for kato-temple, otherwise simple |
[out] | lammin | Lower bound of the spectrum |
[out] | lammax | Upper bound of the spectrum |
[out] | fstats | File stream which stats are printed to |
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().
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().
void randn_double | ( | int | n, |
double * | v | ||
) |
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].
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 |
maxit | max Num of outer Lanczos steps allowed –[max dim of Krylov subspace] |
tol | tolerance 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 |
vinit | initial vector for Lanczos – [optional] |
[out] | rat | A struct containing the rational filter |
[out] | nevOut | Number of eigenvalues/vectors computed |
[out] | Wo | A set of eigenvectors [n x nevOut matrix] |
[out] | lamo | Associated eigenvalues [nevOut x 1 vector] |
[out] | reso | Associated residual norms [nev x 1 vector] |
[out] | fstats | File 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().
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].
[in] | lanm | Dimension of Krylov subspace [restart dimension] |
[in] | nev | Estimate 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] | 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 |
[in] | maxit | max Num of outer Lanczos iterations (restarts) allowed – Each restart may or use the full lanm lanczos steps or fewer. |
[in] | tol | tolerance for convergence. stop when ||res||< tol |
[in] | vinit | initial vector for Lanczos – [optional] |
[in] | rat | a struct containing the parameters of the rational. |
[out] | nev2 | Number of eigenvalues/vectors computed |
[out] | W | A set of eigenvectors [n x nev2 matrix] of unit 2-norm for standard eig prob of unit B-norm for generalized eig prob |
[out] | vals | Associated eigenvalues [nev2 x 1 vector] |
[out] | resW | Associated residual norms [nev x 1 vector] 2-norm for standard eig prob B-norm for generalized eig prob |
[out] | fstats | File stream which stats are printed to |
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().
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().
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().
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().
int SetAMatvec | ( | int | n, |
MVFunc | func, | ||
void * | data | ||
) |
Set the user-input matvec routine and the associated data for A. Save them in evsldata.
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().
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().
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().
int SetBMatvec | ( | int | n, |
MVFunc | func, | ||
void * | data | ||
) |
Set the user-input matvec routine and the associated data for B. Save them in evsldata.
Definition at line 117 of file evsl.c.
References _evsldata::Bmv, Calloc, _EVSLMatvec::data, _EVSLMatvec::func, and _evsldata::n.
Referenced by evsl_setbmv().
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().
void SetDiagScal | ( | double * | 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().
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().
int SetStdEig | ( | ) |
Set the problem to standard eigenvalue problem.
Definition at line 148 of file evsl.c.
References _evsldata::ifGenEv.
Referenced by main().
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().
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().
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().
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).
*sli | see above (output) |
*mu | coeffs of polynomial (input) |
Mdeg | degree of polynomial (input) |
*intv | an 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_int | number of slices wanted (input) |
npts | number 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().
void spslicer2 | ( | double * | xi, |
double * | yi, | ||
int | n_int, | ||
int | npts, | ||
double * | sli | ||
) |
Interval partitioner based for Lanczos DOS output
[in] | xi | coordinates of interval [a b] |
[in] | yi | yi[k] = integral of the does from a to xi[k] |
[in] | n_int | Number of desired sub-intervals |
[in] | npts | number of integration points (length of xi) |
[out] | sli | Array 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().
void StatsPrint | ( | FILE * | fstats | ) |
Definition at line 14 of file stats.c.
References evslstat, _evslstat::n_mvA, _evslstat::n_mvB, _evslstat::n_polAv, _evslstat::n_ratAv, _evslstat::n_svASigB, _evslstat::n_svB, _evslstat::t_blas, _evslstat::t_eig, _evslstat::t_iter, _evslstat::t_mvA, _evslstat::t_mvB, _evslstat::t_polAv, _evslstat::t_ratAv, _evslstat::t_reorth, _evslstat::t_ritz, _evslstat::t_setASigBsv, _evslstat::t_setBsv, _evslstat::t_sth, _evslstat::t_svASigB, and _evslstat::t_svB.
Referenced by main().
void StatsReset | ( | ) |
Definition at line 79 of file stats.c.
References evslstat.
Referenced by EVSLStart(), and main().
int UnsetAMatvec | ( | ) |
int UnsetBMatvec | ( | ) |