EVSL  1.1.0
EigenValues Slicing Library
 All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros
MMPLanR.c
Go to the documentation of this file.
1 #include <stdio.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include <sys/stat.h>
5 #include <math.h>
6 #include "evsl.h"
7 #include "io.h"
8 
9 #define max(a, b) ((a) > (b) ? (a) : (b))
10 #define min(a, b) ((a) < (b) ? (a) : (b))
11 #define TRIV_SLICER 0
12 
13 int main () {
14  int ierr = 0;
15  /*--------------------------------------------------------------
16  * this tests the spectrum slicing idea for a generic matrix
17  * read in sparse matrix format
18  *-------------------------------------------------------------*/
19  int n=0, sl, i, j, mlan, nev, totcnt;
20  //int nnz;
21  double a, b, ecount, xintv[4];
22  double lmin, lmax;
23  double *alleigs;
24  int n_intv; // number of subintervals (slices)
25  int npts; // number of integration points for eigenvalue count
26  //--------------------related to kpmdos
27  //int Mdeg = 180, nvec = 100;
28  /*-------------------- matrix A: coo format and csr format */
29  cooMat Acoo;
30  csrMat Acsr;
31  /* tolerance to accept ev */
32  double tol;
33  /* stats */
34  /* #ev computed; the computed eig val/vec,
35  * residuals of the computed eig pairs */
36  int nevOut;
37  // double *lam=NULL, *Y=NULL, *res=NULL ;
38  /*-------------------- mim/max degree of polynomial, max Lanczos steps */
39  int max_its, Mdeg, nvec;
40  /*-------------------- IO */
41  FILE *flog = stdout, *fmat = NULL, *fstats = NULL;
42  io_t io;
43  int numat, mat;
44  char line[MAX_LINE];
45  /* initial vector: random */
46  double *vinit;
47  polparams pol;
48  //-------------------- tolerance for stopping criterion
49  tol = 1e-8;
50  //-------------------- slicer parameters
51  Mdeg = 300;
52  nvec = 60;
53  /*-------------------- start EVSL */
54  EVSLStart();
55  //-------------------- interior eigensolver parameters
56  double *mu = malloc((Mdeg+1)*sizeof(double));
57  int *counts;
58  double *sli;
59  /*------------------ file "matfile" contains paths to matrices */
60  if ( NULL == ( fmat = fopen( "matfile", "r" ) ) ) {
61  fprintf( flog, "Can't open matfile...\n" );
62  exit(2);
63  }
64  /*-------------------- read number of matrices ..*/
65  memset( line, 0, MAX_LINE );
66  if (NULL == fgets( line, MAX_LINE, fmat )) {
67  fprintf( flog, "error in reading matfile...\n" );
68  exit(2);
69  }
70  if( ( numat = atoi( line ) ) <= 0 ) {
71  fprintf( flog, "Invalid count of matrices...\n" );
72  exit(3);
73  }
74  /*-------------------- LOOP through matrices -*/
75  for(mat = 1; mat <= numat; mat++) {
76  if(get_matrix_info(fmat, &io) != 0) {
77  fprintf(flog, "Invalid format in matfile ...\n");
78  exit(5);
79  }
80  /*----------------input matrix and interval information -*/
81  fprintf(flog, "MATRIX: %s...\n", io.MatNam);
82  a = io.a;
83  b = io.b;
84  n_intv = io.n_intv;
85 
86  struct stat st = {0}; /* Make sure OUT directory exists */
87  if (stat("OUT", &st) == -1) {
88  mkdir("OUT", 0750);
89  }
90 
91  /*-------------------- path to write the output files*/
92  char path[1024];
93  strcpy( path, "OUT/MMPLanR_");
94  strcat( path, io.MatNam);
95  fstats = fopen(path,"w"); // write all the output to the file io.MatNam
96  if (!fstats) {
97  printf(" failed in opening output file in OUT/\n");
98  fstats = stdout;
99  }
100  fprintf(fstats, "MATRIX: %s...\n", io.MatNam);
101  fprintf(fstats,"Partition the interval of interest [%f,%f] into %d slices\n", a,b,n_intv);
102  counts = malloc(n_intv*sizeof(int));
103  sli = malloc((n_intv+1)*sizeof(double));
104  /*-------------------- Read matrix - case: COO/MatrixMarket formats */
105  if (io.Fmt > HB) {
106  ierr =read_coo_MM(io.Fname, 1, 0, &Acoo);
107  if (ierr == 0) {
108  fprintf(fstats,"matrix read successfully\n");
109  //nnz = Acoo.nnz;
110  n = Acoo.nrows;
111  }
112  else {
113  fprintf(flog, "read_coo error = %d\n", ierr);
114  exit(6);
115  }
116  /*-------------------- conversion from COO to CSR format */
117  ierr = cooMat_to_csrMat(0, &Acoo, &Acsr);
118  }
119  if (io.Fmt == HB) {
120  fprintf(flog, "HB FORMAT not supported (yet) * \n");
121  exit(7);
122  }
123  /*-------------------- set the left-hand side matrix A */
124  SetAMatrix(&Acsr);
125  /*-------------------- define ChebLanTr parameters */
126  alleigs = (double *) malloc(n*sizeof(double));
127  vinit = (double *) malloc(n*sizeof(double));
128  rand_double(n, vinit);
129  /*-------------------- get lambda_min lambda_max estimates */
130  ierr = LanTrbounds(50, 200, 1e-8, vinit, 1, &lmin, &lmax, fstats);
131  fprintf(fstats, "Step 0: Eigenvalue bounds for A: [%.15e, %.15e]\n", lmin, lmax);
132  /*-------------------- define [a b] now so we can get estimates now
133  on number of eigenvalues in [a b] from kpmdos */
134  fprintf(fstats," --> interval: a %9.3e b %9.3e \n",a, b);
135  /*-------------------- define kpmdos parameters */
136  xintv[0] = a; xintv[1] = b;
137  xintv[2] = lmin; xintv[3] = lmax;
138  //-------------------- kpmdos or triv_slicer
139  if (TRIV_SLICER) {
140  linspace(a, b, n_intv+1, sli);
141  } else {
142  double t = evsl_timer();
143  ierr = kpmdos(Mdeg, 1, nvec, xintv, mu, &ecount);
144  t = evsl_timer() - t;
145  if (ierr) {
146  printf("kpmdos error %d\n", ierr);
147  return 1;
148  }
149  fprintf(fstats, " Time to build DOS (kpmdos) was : %10.2f \n",t);
150  fprintf(fstats, "Step 1a: Estimated eig count in interval - %.15e \n",ecount);
151  if ((ecount <0) | (ecount > n)) {
152  printf(" e-count estimate is incorrect \n ");
153  exit(7);
154  }
155  /*-------------------- error in ecount */
156  if (ecount < 1 || ecount > n)
157  exit(6);
158  /*-------------------- define slicer parameters */
159  npts = 10 * ecount;
160  ierr = spslicer(sli, mu, Mdeg, xintv, n_intv, npts);
161  if (ierr) {
162  printf("spslicer error %d\n", ierr);
163  return 1;
164  }
165  }
166  fprintf(fstats,"DOS parameters: Mdeg = %d, nvec = %d, npts = %d\n",
167  Mdeg, nvec, npts);
168  fprintf(fstats, "Step 1b: Slices found: \n");
169  for (j=0; j<n_intv;j++) {
170  fprintf(fstats, " %2d: [%.15e , %.15e]\n", j+1, sli[j],sli[j+1]);
171  }
172  //-------------------- # eigs per slice
173  //-------------------- approximate number of eigenvalues wanted
174  double fac = 1.2; // this factor insures that # of eigs per slice is slightly overestimated
175  nev = (int) (1 + ecount / ((double) n_intv)); // # eigs per slice
176  nev = (int)(fac*nev); // want an overestimate of ev_int
177  fprintf(fstats,"Step 2: In each slice compute %d eigenvalues ... \n", nev);
178  /*-------------------- MAIN intv LOOP: for each sclice Do: */
179  totcnt = 0;
180  for (sl=0; sl<n_intv; sl++) {
181  fprintf(fstats,"======================================================\n");
182  double *lam, *Y, *res;
183  int *ind;
184  //--------------------
185  a = sli[sl];
186  b = sli[sl+1];
187  fprintf(fstats, " subinterval: [% 12.4e , % 12.4e]\n", a, b);
188  //-------------------- Parameters for ChebLanTr
189  mlan = max(4*nev,300); mlan = min(mlan, n);
190  max_its = 3*mlan; // max number of Lanczos iterations
191  //fprintf(fstats, "Thick Restarted Lanczos with dimension %d\n", mlan);
192  //fprintf(fstats, "Max Lanczos steps %d\n", max_its);
193  xintv[0] = a;
194  xintv[1] = b;
195  //-------------------- set up default parameters for pol.
196  set_pol_def(&pol);
197  // can change default values here e.g.
198  pol.damping = 2; pol.thresh_int = 0.8;
199  //pol.max_deg = 500;
200  //-------------------- Now determine polymomial
201  find_pol(xintv, &pol);
202  fprintf(fstats, " polynomial [type = %d], deg %d, bar %e gam %e\n",
203  pol.type, pol.deg, pol.bar, pol.gam);
204  //-------------------- Call ChebLanTr
205  ierr = ChebLanTr(mlan, nev, xintv, max_its, tol, vinit,
206  &pol, &nevOut, &lam, &Y, &res, fstats);
207  if (ierr) {
208  printf("ChebLanTr error %d\n", ierr);
209  return 1;
210  }
211  /* sort the eigenvals: ascending order
212  * ind: keep the orginal indices */
213  ind = (int *) malloc(nevOut*sizeof(int));
214  sort_double(nevOut, lam, ind);
215 
216  /* print eigenvalues */
217  fprintf(fstats, "- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -\n");
218  fprintf(fstats, " Eigenvalues in [%f, %f]\n",a,b);
219  fprintf(fstats, "- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -\n");
220  fprintf(fstats, " Computed [%d out of %d estimated] ||Res|| ", nevOut, nev);
221  fprintf(fstats, "\n");
222  for (i=0; i<nevOut; i++) {
223  fprintf(fstats, " % .15e %.1e", lam[i], res[ind[i]]);
224  fprintf(fstats,"\n");
225  }
226  fprintf(fstats, "- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -\n");
227  memcpy(&alleigs[totcnt],lam,nevOut*sizeof(double));
228  totcnt += nevOut;
229  counts[sl] = nevOut;
230  //-------------------- free memory within this loop
231  if (lam) free(lam);
232  if (Y) free(Y);
233  if (res) free(res);
234  free(ind);
235  free_pol(&pol);
236  /*-------------------- end slice loop */
237  }
238  fprintf(fstats," --> Total eigenvalues found = %d\n",totcnt);
239  sprintf(path, "OUT/EigsOut_PLanR_%s", io.MatNam);
240  FILE *fmtout = fopen(path,"w");
241  if (fmtout) {
242  for (j=0; j<totcnt; j++)
243  fprintf(fmtout, "%.15e\n", alleigs[j]);
244  fclose(fmtout);
245  }
246  /*-------------------- free memory */
247  free(counts);
248  free(sli);
249  free(vinit);
250  free_coo(&Acoo);
251  free_csr(&Acsr);
252  free(alleigs);
253  if (fstats != stdout) fclose(fstats);
254  /*-------------------- end matrix loop */
255  }
256  //-------------------- free rest of memory
257  free(mu);
258  if( flog != stdout ) fclose ( flog );
259  fclose( fmat );
260 
261  /*-------------------- finalize EVSL */
262  EVSLFinish();
263 
264  return 0;
265 }
266 
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...
Definition: spslice.c:32
int deg
Definition: struct.h:70
void free_pol(polparams *pol)
Definition: chebpoly.c:488
void free_coo(cooMat *coo)
memory deallocation for coo matrix
Definition: spmat.c:126
double b
Definition: io.h:25
#define MAX_LINE
Definition: io.h:5
#define max(a, b)
Definition: MMPLanR.c:9
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].
Definition: cheblanTr.c:62
int SetAMatrix(csrMat *A)
Set the matrix A.
Definition: evsl.c:68
#define min(a, b)
Definition: MMPLanR.c:10
#define TRIV_SLICER
Definition: MMPLanR.c:11
void rand_double(int n, double *v)
Definition: vect.c:11
int Fmt
Definition: io.h:20
int main()
Definition: MMPLanR.c:18
int cooMat_to_csrMat(int cooidx, cooMat *coo, csrMat *csr)
convert coo to csr
Definition: spmat.c:135
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].
Definition: lanTrbounds.c:38
int find_pol(double *intv, polparams *pol)
Sets the values in pol.
Definition: chebpoly.c:361
double gam
Definition: struct.h:64
void sort_double(int n, double *v, int *ind)
Definition: vect.c:92
int EVSLStart()
Initialize evslData.
Definition: evsl.c:27
void free_csr(csrMat *csr)
memory deallocation for csr matrix
Definition: spmat.c:91
#define HB
Definition: io.h:7
char Fname[MAX_LINE]
Definition: io.h:17
char MatNam[MaxNamLen]
Definition: io.h:18
void set_pol_def(polparams *pol)
set default values for polparams struct.
Definition: chebpoly.c:18
double thresh_int
Definition: struct.h:52
int type
Definition: struct.h:59
sparse matrix format: the compressed sparse row (CSR) format, 0-based
Definition: struct.h:31
This file contains function prototypes and constant definitions for EVSL.
Definition: io.h:14
int EVSLFinish()
Finish EVSL.
Definition: evsl.c:48
int get_matrix_info(FILE *fmat, io_t *pio)
Definition: io.c:22
double a
Definition: io.h:24
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 eac...
Definition: spslice.c:204
int damping
Definition: struct.h:50
int nrows
Definition: struct.h:17
int read_coo_MM(const char *matfile, int idxin, int idxout, cooMat *Acoo)
Definition: io.c:66
parameters for polynomial filter
Definition: struct.h:45
void linspace(double a, double b, int num, double *arr)
Definition: vect.c:54
sparse matrix format: the coordinate (COO) format, 0-based
Definition: struct.h:16
double evsl_timer()
evsl timer for mac
Definition: mactime.c:14
double bar
Definition: struct.h:65
int n_intv
Definition: io.h:23