EVSL  1.1.0
EigenValues Slicing Library
 All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros
LapRLanN.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 <complex.h>
7 #include "evsl.h"
8 #include "io.h"
9 #include "evsl_direct.h"
10 
11 #define max(a, b) ((a) > (b) ? (a) : (b))
12 #define min(a, b) ((a) < (b) ? (a) : (b))
13 
14 int findarg(const char *argname, ARG_TYPE type, void *val, int argc, char **argv);
15 int lapgen(int nx, int ny, int nz, cooMat *Acoo);
16 int exeiglap3(int nx, int ny, int nz, double a, double b, int *m, double **vo);
17 
18 int main(int argc, char *argv[]) {
19  /*------------------------------------------------------------
20  generates a laplacean matrix on an nx x ny x nz mesh
21  and computes all eigenvalues in a given interval [a b]
22  The default set values are
23  nx = 41; ny = 53; nz = 1;
24  a = 0.4; b = 0.8;
25  nslices = 1 [one slice only]
26  other parameters
27  tol [tolerance for stopping - based on residual]
28  Mdeg = pol. degree used for DOS
29  nvec = number of sample vectors used for DOS
30  This uses:
31  Non-restart Lanczos with rational filtering
32  ------------------------------------------------------------*/
33  int n, nx, ny, nz, i, j, npts, nslices, nvec, Mdeg, nev,
34  max_its, ev_int, sl, flg, ierr;
35  /* find the eigenvalues of A in the interval [a,b] */
36  double a, b, lmax, lmin, ecount, tol, *sli, *mu;
37  double xintv[4];
38  /* initial vector: random */
39  double *vinit;
40  /* parameters for rational filter */
41  int num = 1; // number of poles used for each slice
42  int pow = 2; // multiplicity of each pole
43  double beta = 0.01; // beta in the LS approximation
44 
45  struct stat st = {0}; /* Make sure OUT directory exists */
46  if (stat("OUT", &st) == -1) {
47  mkdir("OUT", 0750);
48  }
49 
50  FILE *fstats = NULL;
51  if (!(fstats = fopen("OUT/LapRLanN.out","w"))) {
52  printf(" failed in opening output file in OUT/\n");
53  fstats = stdout;
54  }
55  /*-------------------- matrix A: coo format and csr format */
56  cooMat Acoo;
57  csrMat Acsr;
58 #if CXSPARSE == 1
59  printf("-----------------------------------------\n");
60  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");
61  printf("-----------------------------------------\n");
62 #endif
63  /*-------------------- default values */
64  nx = 41;
65  ny = 53;
66  nz = 1;
67  a = 0.4;
68  b = 0.8;
69  nslices = 4;
70  //-----------------------------------------------------------------------
71  //-------------------- reset some default values from command line [Yuanzhe/]
72  /* user input from command line */
73  flg = findarg("help", NA, NULL, argc, argv);
74  if (flg) {
75  printf("Usage: ./testL.ex -nx [int] -ny [int] -nz [int] -a [double] -b [double] -nslices [int]\n");
76  return 0;
77  }
78  findarg("nx", INT, &nx, argc, argv);
79  findarg("ny", INT, &ny, argc, argv);
80  findarg("nz", INT, &nz, argc, argv);
81  findarg("a", DOUBLE, &a, argc, argv);
82  findarg("b", DOUBLE, &b, argc, argv);
83  findarg("nslices", INT, &nslices, argc, argv);
84  fprintf(fstats,"used nx = %3d ny = %3d nz = %3d",nx,ny,nz);
85  fprintf(fstats," [a = %4.2f b= %4.2f], nslices=%2d \n",a,b,nslices);
86  //-------------------- eigenvalue bounds set by hand.
87  lmin = 0.0;
88  lmax = ((nz == 1)? 8.0 :12) ;
89  xintv[0] = a;
90  xintv[1] = b;
91  xintv[2] = lmin;
92  xintv[3] = lmax;
93  tol = 1e-6;
94  n = nx * ny * nz;
95  /*-------------------- generate 2D/3D Laplacian matrix
96  * saved in coo format */
97  ierr = lapgen(nx, ny, nz, &Acoo);
98  /*-------------------- convert coo to csr */
99  ierr = cooMat_to_csrMat(0, &Acoo, &Acsr);
100  /* output the problem settings */
101  fprintf(fstats, "- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -\n");
102  fprintf(fstats, "Laplacian: %d x %d x %d, n = %d, nnz = %d\n", nx, ny, nz, n, Acoo.nnz);
103  fprintf(fstats, "Interval: [%20.15f, %20.15f] -- %d slices \n", a, b, nslices);
104  fprintf(fstats, "Step 0: Eigenvalue bound s for A: [%.15e, %.15e]\n", lmin, lmax);
105  /*-------------------- call kpmdos to get the DOS for dividing the spectrum*/
106  /*-------------------- define kpmdos parameters */
107  Mdeg = 300;
108  nvec = 60;
109  /*-------------------- start EVSL */
110  EVSLStart();
111  /*-------------------- set the left-hand side matrix A */
112  SetAMatrix(&Acsr);
113  //-------------------- call kpmdos
114  mu = malloc((Mdeg+1)*sizeof(double));
115  double t = evsl_timer();
116  ierr = kpmdos(Mdeg, 1, nvec, xintv, mu, &ecount);
117  t = evsl_timer() - t;
118  if (ierr) {
119  printf("kpmdos error %d\n", ierr);
120  return 1;
121  }
122  fprintf(fstats, " Time to build DOS (kpmdos) was : %10.2f \n",t);
123  fprintf(fstats, " estimated eig count in interval: %10.2e \n",ecount);
124  //-------------------- call splicer to slice the spectrum
125  npts = 10 * ecount;
126  sli = malloc((nslices+1)*sizeof(double));
127  ierr = spslicer(sli, mu, Mdeg, xintv, nslices, npts);
128  if (ierr) {
129  printf("spslicer error %d\n", ierr);
130  return 1;
131  }
132  printf("==================== SLICES FOUND ====================\n");
133  for (j=0; j<nslices;j++) {
134  printf(" %2d: [% .15e , % .15e]\n", j+1, sli[j],sli[j+1]);
135  }
136  //-------------------- # eigs per slice
137  ev_int = (int) (1 + ecount / ((double) nslices));
138  //-------------------- initial vector
139  vinit = (double*) malloc(n*sizeof(double));
140  rand_double(n, vinit);
141  //-------------------- For each slice call RatLanrNr
142  for (sl=0; sl<nslices; sl++) {
143  printf("======================================================\n");
144  int nev2;
145  double *lam, *Y, *res;
146  int *ind;
147  int nev_ex;
148  double *lam_ex;
149  //--------------------
150  StatsReset();
151  a = sli[sl];
152  b = sli[sl+1];
153  printf(" subinterval: [%.4e , %.4e]\n", a, b);
154  double intv[4];
155  intv[0] = a;
156  intv[1] = b;
157  intv[2] = lmin;
158  intv[3] = lmax;
159  // find the rational filter on this slice
160  ratparams rat;
161  // setup default parameters for rat
162  set_ratf_def(&rat);
163  // change some default parameters here:
164  rat.pw = pow;
165  rat.num = num;
166  rat.beta = beta;
167  // now determine rational filter
168  find_ratf(intv, &rat);
169  /*------------ use direct solver function */
170  void **solshiftdata = (void **) malloc(rat.num*sizeof(void *));
171  /*------------ factoring the shifted matrices and store the factors */
172  SetupASIGMABSolDirect(&Acsr, NULL, rat.num, rat.zk, solshiftdata);
173  /*------------ give the data to rat */
174  SetASigmaBSol(&rat, NULL, ASIGMABSolDirect, solshiftdata);
175  //-------------------- approximate number of eigenvalues wanted
176  nev = ev_int+2;
177  //-------------------- maximal Lanczos iterations
178  max_its = max(4*nev,300); max_its = min(max_its, n);
179  //-------------------- RationalLanNr
180  ierr = RatLanNr(intv, max_its, tol, vinit, &rat, &nev2, &lam,
181  &Y, &res, fstats);
182  if (ierr) {
183  printf("RatLanNr error %d\n", ierr);
184  return 1;
185  }
186  /* [compute residual] already computed in res */
187  /* sort the eigenvals: ascending order
188  * ind: keep the orginal indices */
189  ind = (int *) malloc(nev2*sizeof(int));
190  sort_double(nev2, lam, ind);
191  /* compute exact eigenvalues */
192  exeiglap3(nx, ny, nz, a, b, &nev_ex, &lam_ex);
193  printf(" number of eigenvalues: %d, found: %d\n", nev_ex, nev2);
194  /* print eigenvalues */
195  fprintf(fstats, " Eigenvalues in [a, b]\n");
196  fprintf(fstats, " Computed [%d] ||Res|| Exact [%d]",
197  nev2, nev_ex);
198  if (nev2 == nev_ex) {
199  fprintf(fstats, " Err");
200  }
201  fprintf(fstats, "\n");
202  for (i=0; i<max(nev2, nev_ex); i++) {
203  if (i < nev2) {
204  fprintf(fstats, "% .15e %.1e", lam[i], res[ind[i]]);
205  } else {
206  fprintf(fstats, " ");
207  }
208  if (i < nev_ex) {
209  fprintf(fstats, " % .15e", lam_ex[i]);
210  }
211  if (nev2 == nev_ex) {
212  fprintf(fstats, " % .1e", lam[i]-lam_ex[i]);
213  }
214  fprintf(fstats,"\n");
215  if (i>50) {
216  fprintf(fstats," -- More not shown --\n");
217  break;
218  }
219  }
220  //-------------------- free allocated space withing this scope
221  if (lam) free(lam);
222  if (Y) free(Y);
223  if (res) free(res);
224  free(ind);
225  free(lam_ex);
226  FreeASIGMABSolDirect(rat.num, solshiftdata);
227  free(solshiftdata);
228  free_rat(&rat);
229  StatsPrint(fstats);
230  } //for (sl=0; sl<nslices; sl++)
231  //-------------------- free other allocated space
232  free(vinit);
233  free(sli);
234  free_coo(&Acoo);
235  free_csr(&Acsr);
236  free(mu);
237  fclose(fstats);
238  /*-------------------- finalize EVSL */
239  EVSLFinish();
240  return 0;
241 }
242 
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
void free_coo(cooMat *coo)
memory deallocation for coo matrix
Definition: spmat.c:126
void ASIGMABSolDirect(int n, double *br, double *bi, double *xr, double *xz, void *data)
complex linear solver routine passed to evsl
complex double * zk
Definition: struct.h:126
int num
Definition: struct.h:114
void free_rat(ratparams *rat)
Definition: ratfilter.c:424
ARG_TYPE
Definition: io.h:33
int SetupASIGMABSolDirect(csrMat *A, csrMat *BB, int num, complex double *zk, void **data)
setup CXsparse solver for A - SIGMA B
int nnz
Definition: struct.h:17
Definition: io.h:35
int SetAMatrix(csrMat *A)
Set the matrix A.
Definition: evsl.c:68
void StatsPrint(FILE *fstats)
Definition: stats.c:14
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
void FreeASIGMABSolDirect(int num, void **data)
free the data needed by CXSparse
int main(int argc, char *argv[])
Definition: LapRLanN.c:18
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
void sort_double(int n, double *v, int *ind)
Definition: vect.c:92
int EVSLStart()
Initialize evslData.
Definition: evsl.c:27
int pw
Definition: struct.h:115
void free_csr(csrMat *csr)
memory deallocation for csr matrix
Definition: spmat.c:91
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].
Definition: ratlanNr.c:54
int lapgen(int nx, int ny, int nz, cooMat *Acoo)
Laplacean Matrix generator.
Definition: lapl.c:17
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 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
void StatsReset()
Definition: stats.c:79
double beta
Definition: struct.h:117
void set_ratf_def(ratparams *rat)
Sets default values for ratparams struct.
Definition: ratfilter.c:345
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 ...
Definition: evsl.c:168
Definitions used for direct solver interface.
int findarg(const char *argname, ARG_TYPE type, void *val, int argc, char **argv)
Definition: io.c:171
int find_ratf(double *intv, ratparams *rat)
Definition: ratfilter.c:375
sparse matrix format: the coordinate (COO) format, 0-based
Definition: struct.h:16
#define max(a, b)
Definition: LapRLanN.c:11
double evsl_timer()
evsl timer for mac
Definition: mactime.c:14
parameters for rational filter
Definition: struct.h:112
Definition: io.h:37
#define min(a, b)
Definition: LapRLanN.c:12