EVSL  1.1.0
EigenValues Slicing Library
 All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros
LapPLanR.c
Go to the documentation of this file.
1 #include <stdlib.h>
2 #include <stdio.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 int findarg(const char *argname, ARG_TYPE type, void *val, int argc, char **argv);
12 int lapgen(int nx, int ny, int nz, cooMat *Acoo);
13 int exeiglap3(int nx, int ny, int nz, double a, double b, int *m, double **vo);
14 
15 int main(int argc, char *argv[]) {
16  /*------------------------------------------------------------
17  generates a laplacean matrix on an nx x ny x nz mesh
18  and computes all eigenvalues in a given interval [a b]
19  The default set values are
20  nx = 41; ny = 53; nz = 1;
21  a = 0.4; b = 0.8;
22  nslices = 1 [one slice only]
23  other parameters
24  tol [tolerance for stopping - based on residual]
25  Mdeg = pol. degree used for DOS
26  nvec = number of sample vectors used for DOS
27  This uses:
28  Thick-restart Lanczos with polynomial filtering
29  ------------------------------------------------------------*/
30  int n, nx, ny, nz, i, j, npts, nslices, nvec, Mdeg, nev,
31  mlan, max_its, ev_int, sl, flg, ierr;
32  /* find the eigenvalues of A in the interval [a,b] */
33  double a, b, lmax, lmin, ecount, tol, *sli, *mu;
34  double xintv[4];
35  double *vinit;
36  polparams pol;
37 
38  struct stat st = {0}; /* Make sure OUT directory exists */
39  if (stat("OUT", &st) == -1) {
40  mkdir("OUT", 0750);
41  }
42 
43  FILE *fstats = NULL;
44  if (!(fstats = fopen("OUT/LapPLanR.out","w"))) {
45  printf(" failed in opening output file in OUT/\n");
46  fstats = stdout;
47  }
48  /*-------------------- matrix A: coo format and csr format */
49  cooMat Acoo;
50  csrMat Acsr;
51  /*-------------------- default values */
52  nx = 41;
53  ny = 53;
54  nz = 8;
55  a = 0.4;
56  b = 0.8;
57  nslices = 4;
58  //-----------------------------------------------------------------------
59  //-------------------- reset some default values from command line
60  /* user input from command line */
61  flg = findarg("help", NA, NULL, argc, argv);
62  if (flg) {
63  printf("Usage: ./testL.ex -nx [int] -ny [int] -nz [int] -a [double] -b [double] -nslices [int]\n");
64  return 0;
65  }
66  findarg("nx", INT, &nx, argc, argv);
67  findarg("ny", INT, &ny, argc, argv);
68  findarg("nz", INT, &nz, argc, argv);
69  findarg("a", DOUBLE, &a, argc, argv);
70  findarg("b", DOUBLE, &b, argc, argv);
71  findarg("nslices", INT, &nslices, argc, argv);
72  fprintf(fstats,"used nx = %3d ny = %3d nz = %3d",nx,ny,nz);
73  fprintf(fstats," [a = %4.2f b= %4.2f], nslices=%2d \n",a,b,nslices);
74  //-------------------- eigenvalue bounds set by hand.
75  lmin = 0.0;
76  lmax = nz == 1 ? 8.0 : 12.0;
77  xintv[0] = a;
78  xintv[1] = b;
79  xintv[2] = lmin;
80  xintv[3] = lmax;
81  tol = 1e-8;
82  n = nx * ny * nz;
83  /*-------------------- output the problem settings */
84  fprintf(fstats, "- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -\n");
85  fprintf(fstats, "Laplacian: %d x %d x %d, n = %d\n", nx, ny, nz, n);
86  fprintf(fstats, "Interval: [%20.15f, %20.15f] -- %d slices \n", a, b, nslices);
87  /*-------------------- generate 2D/3D Laplacian matrix
88  * saved in coo format */
89  ierr = lapgen(nx, ny, nz, &Acoo);
90  /*-------------------- convert coo to csr */
91  ierr = cooMat_to_csrMat(0, &Acoo, &Acsr);
92  /*-------------------- step 0: get eigenvalue bounds */
93  fprintf(fstats, "Step 0: Eigenvalue bound s for A: [%.15e, %.15e]\n", lmin, lmax);
94  /*-------------------- call kpmdos to get the DOS for dividing the spectrum*/
95  /*-------------------- define kpmdos parameters */
96  Mdeg = 300;
97  nvec = 60;
98  /*-------------------- start EVSL */
99  EVSLStart();
100  /*-------------------- set the left-hand side matrix A */
101  SetAMatrix(&Acsr);
102  /*-------------------- call kpmdos */
103  mu = malloc((Mdeg+1)*sizeof(double));
104  double t = evsl_timer();
105  ierr = kpmdos(Mdeg, 1, nvec, xintv, mu, &ecount);
106  t = evsl_timer() - t;
107  if (ierr) {
108  printf("kpmdos error %d\n", ierr);
109  return 1;
110  }
111  fprintf(fstats, " Time to build DOS (kpmdos) was : %10.2f \n",t);
112  fprintf(fstats, " estimated eig count in interval: %.15e \n",ecount);
113  //-------------------- call splicer to slice the spectrum
114  npts = 10 * ecount;
115  sli = malloc((nslices+1)*sizeof(double));
116 
117  fprintf(fstats,"DOS parameters: Mdeg = %d, nvec = %d, npnts = %d\n",Mdeg, nvec, npts);
118  ierr = spslicer(sli, mu, Mdeg, xintv, nslices, npts);
119  if (ierr) {
120  printf("spslicer error %d\n", ierr);
121  return 1;
122  }
123  printf("==================== SLICES FOUND ====================\n");
124  for (j=0; j<nslices;j++) {
125  printf(" %2d: [% .15e , % .15e]\n", j+1, sli[j],sli[j+1]);
126  }
127  //-------------------- # eigs per slice
128  ev_int = (int) (1 + ecount / ((double) nslices));
129  //-------------------- initial vector
130  vinit = (double *) malloc(n*sizeof(double));
131  rand_double(n, vinit);
132  //-------------------- debug only :
133  // save_vec(n, vinit, "OUT/vinit.mtx");
134  //-------------------- For each slice call ChebLanr
135  for (sl =0; sl<nslices; sl++){
136  printf("======================================================\n");
137  int nev2;
138  double *lam, *Y, *res;
139  int *ind;
140  int nev_ex;
141  double *lam_ex;
142  //--------------------
143  a = sli[sl];
144  b = sli[sl+1];
145  printf(" subinterval: [%.4e , %.4e]\n", a, b);
146  //-------------------- approximate number of eigenvalues wanted
147  nev = ev_int+2;
148  //-------------------- Dimension of Krylov subspace
149  mlan = max(4*nev, 300);
150  mlan = min(mlan, n);
151  max_its = 3*mlan;
152  //-------------------- ChebLanTr
153  xintv[0] = a; xintv[1] = b;
154  xintv[2] = lmin; xintv[3] = lmax;
155  //-------------------- set up default parameters for pol.
156  set_pol_def(&pol);
157  //-------------------- this is to show how you can reset some of the
158  // parameters to determine the filter polynomial
159  pol.damping = 0;
160  //-------------------- use a stricter requirement for polynomial
161  pol.thresh_int = 0.8;
162  pol.thresh_ext = 0.2;
163  pol.max_deg = 10000;
164  // pol.deg = 20 //<< this will force this exact degree . not recommended
165  // it is better to change the values of the thresholds
166  // pol.thresh_ext and plot.thresh_int
167  //-------------------- Now determine polymomial to use
168  find_pol(xintv, &pol);
169 
170  fprintf(fstats, " polynomial [type = %d], deg %d, bar %e gam %e\n",
171  pol.type, pol.deg, pol.bar, pol.gam);
172  //-------------------- then call ChenLanTr
173  ierr = ChebLanTr(mlan, nev, xintv, max_its, tol, vinit,
174  &pol, &nev2, &lam, &Y, &res, fstats);
175  if (ierr) {
176  printf("ChebLanTr error %d\n", ierr);
177  return 1;
178  }
179  /*--------------------- residuals were already computed in res */
180  /* sort the eigenvals: ascending order
181  * ind: keep the orginal indices */
182  ind = (int *) malloc(nev2*sizeof(int));
183  sort_double(nev2, lam, ind);
184  /* compute exact eigenvalues */
185  exeiglap3(nx, ny, nz, a, b, &nev_ex, &lam_ex);
186  printf(" number of eigenvalues: %d, found: %d\n", nev_ex, nev2);
187 
188  /* print eigenvalues */
189  fprintf(fstats, " Eigenvalues in [a, b]\n");
190  fprintf(fstats, " Computed [%d] ||Res|| Exact [%d]", nev2, nev_ex);
191  if (nev2 == nev_ex) {
192  fprintf(fstats, " Err");
193  }
194  fprintf(fstats, "\n");
195  for (i=0; i<max(nev2, nev_ex); i++) {
196  if (i < nev2) {
197  fprintf(fstats, "% .15e %.1e", lam[i], res[ind[i]]);
198  } else {
199  fprintf(fstats, " ");
200  }
201  if (i < nev_ex) {
202  fprintf(fstats, " % .15e", lam_ex[i]);
203  }
204  if (nev2 == nev_ex) {
205  fprintf(fstats, " % .1e", lam[i]-lam_ex[i]);
206  }
207  fprintf(fstats,"\n");
208  if (i>50) {
209  fprintf(fstats," -- More not shown --\n");
210  break;
211  }
212  }
213  //-------------------- free allocated space withing this scope
214  if (lam) free(lam);
215  if (Y) free(Y);
216  if (res) free(res);
217  free_pol(&pol);
218  free(ind);
219  free(lam_ex);
220  }
221  //-------------------- free other allocated space
222  free(vinit);
223  free(sli);
224  free_coo(&Acoo);
225  free_csr(&Acsr);
226  free(mu);
227  fclose(fstats);
228  /*-------------------- finalize EVSL */
229  EVSLFinish();
230  return 0;
231 }
232 
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
ARG_TYPE
Definition: io.h:33
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
Definition: io.h:35
int SetAMatrix(csrMat *A)
Set the matrix A.
Definition: evsl.c:68
void rand_double(int n, double *v)
Definition: vect.c:11
int cooMat_to_csrMat(int cooidx, cooMat *coo, csrMat *csr)
convert coo to csr
Definition: spmat.c:135
Definition: io.h:34
int lapgen(int nx, int ny, int nz, cooMat *Acoo)
Laplacean Matrix generator.
Definition: lapl.c:17
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
int exeiglap3(int nx, int ny, int nz, double a, double b, int *m, double **vo)
Exact eigenvalues of Laplacean in interval [a b].
Definition: lapl.c:83
double thresh_ext
Definition: struct.h:51
void set_pol_def(polparams *pol)
set default values for polparams struct.
Definition: chebpoly.c:18
double thresh_int
Definition: struct.h:52
int findarg(const char *argname, ARG_TYPE type, void *val, int argc, char **argv)
Definition: io.c:171
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.
int max_deg
Definition: struct.h:48
int EVSLFinish()
Finish EVSL.
Definition: evsl.c:48
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
parameters for polynomial filter
Definition: struct.h:45
#define max(a, b)
Definition: LapPLanR.c:9
#define min(a, b)
Definition: LapPLanR.c:10
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 main(int argc, char *argv[])
Definition: LapPLanR.c:15
Definition: io.h:37