10 #define max(a, b) ((a) > (b) ? (a) : (b))
11 #define min(a, b) ((a) < (b) ? (a) : (b))
24 int n = 0, i, j, npts, nslices, nvec, nev, mlan, max_its, ev_int, sl, ierr, totcnt;
26 double a, b, lmax, lmin, ecount, tol, *sli;
40 FILE *flog = stdout, *fmat = NULL;
48 printf(
"-----------------------------------------\n");
49 printf(
"Note: You are using CXSparse for the direct solver. \n We recommend a more performance based direct solver for anything more than basic tests. \n SuiteSparse is supported with a makefile change. \n Using SuiteSparse can result in magnitudes faster times. \n\n");
50 printf(
"-----------------------------------------\n");
57 if (NULL == (fmat = fopen(
"matfile",
"r"))) {
58 fprintf(flog,
"Can't open matfile...\n");
63 if (NULL == fgets(line,
MAX_LINE, fmat)) {
64 fprintf(flog,
"error in reading matfile...\n");
67 if ((numat = atoi(line)) <= 0) {
68 fprintf(flog,
"Invalid count of matrices...\n");
71 for (mat = 1; mat <= numat; mat++) {
73 fprintf(flog,
"Invalid format in matfile ...\n");
77 fprintf(flog,
"MATRIX A: %s...\n", io.
MatNam1);
78 fprintf(flog,
"MATRIX B: %s...\n", io.
MatNam2);
84 if (stat(
"OUT", &st) == -1) {
89 strcpy(path,
"OUT/Lan_MMPLanR_");
91 fstats = fopen(path,
"w");
93 printf(
" failed in opening output file in OUT/\n");
96 fprintf(fstats,
"MATRIX A: %s...\n", io.
MatNam1);
97 fprintf(fstats,
"MATRIX B: %s...\n", io.
MatNam2);
98 fprintf(fstats,
"Partition the interval of interest [%f,%f] into %d slices\n",
100 counts = malloc(nslices *
sizeof(
int));
101 sli = malloc((nslices + 1) *
sizeof(
double));
106 fprintf(fstats,
"matrix read successfully\n");
112 fprintf(flog,
"read_coo error for A = %d\n", ierr);
117 fprintf(fstats,
"matrix read successfully\n");
118 if (Bcoo.
nrows != n) {
122 fprintf(flog,
"read_coo error for B = %d\n", ierr);
128 }
else if (io.
Fmt ==
HB) {
129 fprintf(flog,
"HB FORMAT not supported (yet) * \n");
132 alleigs = malloc(n *
sizeof(
double));
146 vinit = (
double *)malloc(n *
sizeof(
double));
148 ierr =
LanTrbounds(50, 200, 1e-12, vinit, 1, &lmin, &lmax, fstats);
149 fprintf(fstats,
"Step 0: Eigenvalue bound s for B^{-1}*A: [%.15e, %.15e]\n",
159 double *xdos = (
double *)calloc(npts,
sizeof(
double));
160 double *ydos = (
double *)calloc(npts,
sizeof(
double));
161 ierr =
LanDosG(nvec, msteps, npts, xdos, ydos, &ecount, xintv);
164 printf(
"Landos error %d\n", ierr);
167 fprintf(fstats,
" Time to build DOS (Landos) was : %10.2f \n", t);
168 fprintf(fstats,
" estimated eig count in interval: %.15e \n", ecount);
170 fprintf(fstats,
"DOS parameters: msteps = %d, nvec = %d, npnts = %d\n",
172 spslicer2(xdos, ydos, nslices, npts, sli);
173 printf(
"==================== SLICES FOUND ====================\n");
174 for (j = 0; j < nslices; j++) {
175 printf(
" %2d: [% .15e , % .15e]\n", j + 1, sli[j], sli[j + 1]);
178 ev_int = (int)(1 + ecount / ((
double)nslices));
181 for (sl = 0; sl < nslices; sl++) {
182 printf(
"======================================================\n");
184 double *lam, *Y, *res;
189 printf(
" subinterval: [%.15e , %.15e]\n", a, b);
203 fprintf(fstats,
" polynomial deg %d, bar %.15e gam %.15e\n", pol.
deg,
209 mlan =
max(4 * nev, 100);
213 ierr =
ChebLanTr(mlan, nev, xintv, max_its, tol, vinit, &pol, &nev2, &lam,
216 printf(
"ChebLanTr error %d\n", ierr);
222 ind = (
int *)malloc(nev2 *
sizeof(
int));
224 printf(
" number of eigenvalues found: %d\n", nev2);
226 fprintf(fstats,
" Eigenvalues in [a, b]\n");
227 fprintf(fstats,
" Computed [%d] ||Res||\n", nev2);
228 for (i = 0; i < nev2; i++) {
229 fprintf(fstats,
"% .15e %.1e\n", lam[i], res[ind[i]]);
231 fprintf(fstats,
"- - - - - - - - - - - - - - - - - - - - - - - - - - - - "
232 "- - - - - - - - - - - - - - - - - -\n");
233 memcpy(&alleigs[totcnt], lam, nev2 *
sizeof(
double));
247 fprintf(fstats,
" --> Total eigenvalues found = %d\n", totcnt);
248 sprintf(path,
"OUT/EigsOut_Lan_MMPLanR_(%s_%s)", io.
MatNam1, io.
MatNam2);
249 FILE *fmtout = fopen(path,
"w");
251 for (j = 0; j < totcnt; j++)
252 fprintf(fmtout,
"%.15e\n", alleigs[j]);
266 if (fstats != stdout) {
271 if (flog != stdout) {
void free_pol(polparams *pol)
void free_coo(cooMat *coo)
memory deallocation for coo matrix
int SetGenEig()
Set the problem to generalized eigenvalue problem.
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].
int SetAMatrix(csrMat *A)
Set the matrix A.
int SetupBSolDirect(csrMat *B, void **data)
Setup the B-sol by computing the Cholesky factorization of B.
void rand_double(int n, double *v)
int cooMat_to_csrMat(int cooidx, cooMat *coo, csrMat *csr)
convert coo to csr
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].
int find_pol(double *intv, polparams *pol)
Sets the values in pol.
int SetBMatrix(csrMat *B)
Set the B matrix.
int read_coo_MM(const char *matfile, int idxin, int idxout, cooMat *Acoo)
void sort_double(int n, double *v, int *ind)
void spslicer2(double *xi, double *yi, int n_int, int npts, double *sli)
int EVSLStart()
Initialize evslData.
void FreeBSolDirectData(void *data)
Free solver data.
int get_matrix_info(FILE *fmat, io_t *pio)
void free_csr(csrMat *csr)
memory deallocation for csr matrix
void set_pol_def(polparams *pol)
set default values for polparams struct.
int LanDosG(const int nvec, const int msteps, int npts, double *xdos, double *ydos, double *neig, const double *const intv)
int SetLTSol(SolFuncR func, void *data)
Set the solve routine for LT.
sparse matrix format: the compressed sparse row (CSR) format, 0-based
This file contains function prototypes and constant definitions for EVSL.
int EVSLFinish()
Finish EVSL.
void LTSolDirect(double *b, double *x, void *data)
Solver function of L^{T} x = L^{-T}*b.
int SetBSol(SolFuncR func, void *data)
Set the solve routine and the associated data for B.
parameters for polynomial filter
void BSolDirect(double *b, double *x, void *data)
Solver function of B.
Definitions used for direct solver interface.
sparse matrix format: the coordinate (COO) format, 0-based
double evsl_timer()
evsl timer for mac