19 #define FILTER_VINIT 1
61 int RatLanTr(
int lanm,
int nev,
double *intv,
int maxit,
62 double tol,
double *vinit,
ratparams *rat,
int *nev2,
63 double **vals,
double **W,
double **resW, FILE *fstats) {
81 size_t lanm1_l = lanm1;
84 maxit =
min(maxit, n);
96 double done=1.0, dmone=-1.0, dzero=0.0;
98 int Ntest =
min(lanm, nev+50);
101 int i, ll, jl, last_jl;
105 if (check_intv(intv, fstats) < 0) {
107 *vals = NULL; *W = NULL; *resW = NULL;
118 fprintf(fstats,
" Rational-LanTR, dim %d cycle %d \n",lanm, cycle);
123 int min_inner_step = 5;
128 Malloc(V, n_l*lanm1,
double);
132 Malloc(Z, n_l*lanm1,
double);
137 Calloc(T, lanm1_l*lanm1_l,
double);
140 double *Y, *Q, *Lam, *res;
141 Malloc(Y, n_l*nev,
double);
146 Malloc(Q, n_l*nev,
double);
153 int trlen = 0, prtrlen=-1;
155 double *Rval, *Rvec, *resi, *Bvec = NULL;
156 Malloc(Rval, lanm,
double);
157 Malloc(resi, lanm,
double);
158 Malloc(Rvec, n_l*lanm,
double);
160 Malloc(Bvec, n_l*lanm,
double);
166 Malloc(EvecT, lanm1_l*lanm1_l,
double);
171 double *work, *vrand = NULL;
172 size_t work_size = ifGenEv ? 6*n_l : 4*n_l;
173 Malloc(work, work_size,
double);
179 DCOPY(&n, V, &one, Z, &one);
183 DCOPY(&n, vinit, &one, Z, &one);
191 t = 1.0 / sqrt(
DDOT(&n, V, &one, Z, &one));
192 DSCAL(&n, &t, Z, &one);
195 t = 1.0 /
DNRM2(&n, V, &one);
198 DSCAL(&n, &t, V, &one);
205 double beta = 0.0, r, res0;
207 int k = trlen > 0 ? trlen + 1 : 0;
210 if (k+min_inner_step > lanm1) {
212 fprintf(fstats,
"Krylov dim too small for this problem. Try a larger dim\n");
220 double *v = &V[k1*n_l];
221 double *z = &Z[k1*n_l];
223 double *vnew = v + n;
224 double *znew = z + n;
239 for (i=0; i<trlen; i++) {
240 T[i*lanm1_l+i] = Rval[i];
244 s[k1] =
DDOT(&n, v, &one, znew, &one);
246 DGEMV(&cN, &n, &k, &dmone, Z, &n, s, &one, &done, znew, &one);
249 for (i=0; i<k1; i++) {
250 T[trlen*lanm1_l+i] = s[i];
251 T[i*lanm1_l+trlen] = s[i];
252 wn += 2.0 * fabs(s[i]);
254 T[trlen*lanm1_l+trlen] = s[k1];
258 matvec_B(znew, vnew);
260 beta = sqrt(
DDOT(&n, vnew, &one, znew, &one));
263 beta =
DNRM2(&n, vnew, &one);
270 fprintf(fstats,
"it %4d: Lucky breakdown, beta = %.15e\n", k, beta);
285 matvec_B(znew, vnew);
286 beta = sqrt(
DDOT(&n, vnew, &one, znew, &one));
287 double ibeta = 1.0 / beta;
288 DSCAL(&n, &ibeta, vnew, &one);
289 DSCAL(&n, &ibeta, znew, &one);
297 double ibeta = 1.0 / beta;
298 DSCAL(&n, &ibeta, vnew, &one);
303 double ibeta = 1.0 / beta;
304 DSCAL(&n, &ibeta, vnew, &one);
306 DSCAL(&n, &ibeta, znew, &one);
310 T[k1*lanm1_l+k] = beta;
311 T[k*lanm1_l+k1] = beta;
315 Ntest =
max(20,nev-lock+10);
316 last_jl = 0; last_tr = 0.0;
320 double *zold = k > 0 ? Z+(k-1)*n_l : NULL;
324 while (k < lanm && it < maxit) {
327 double *v = &V[(k-1)*n_l];
328 double *z = &Z[(k-1)*n_l];
330 double *vnew = v + n;
331 double *znew = z + n;
347 double nbeta = -beta;
348 DAXPY(&n, &nbeta, zold, &one, znew, &one);
351 double alpha =
DDOT(&n, v, &one, znew, &one);
353 T[(k-1)*lanm1_l+(k-1)] = alpha;
356 double nalpha = -alpha;
357 DAXPY(&n, &nalpha, z, &one, znew, &one);
363 matvec_B(znew, vnew);
365 beta = sqrt(
DDOT(&n, vnew, &one, znew, &one));
378 fprintf(fstats,
"it %4d: Lucky breakdown, beta = %.15e\n", it, beta);
393 matvec_B(znew, vnew);
394 beta = sqrt(
DDOT(&n, vnew, &one, znew, &one));
395 double ibeta = 1.0 / beta;
396 DSCAL(&n, &ibeta, vnew, &one);
397 DSCAL(&n, &ibeta, znew, &one);
405 double ibeta = 1.0 / beta;
406 DSCAL(&n, &ibeta, vnew, &one);
411 double ibeta = 1.0 / beta;
412 DSCAL(&n, &ibeta, vnew, &one);
415 DSCAL(&n, &ibeta, znew, &one);
419 T[k*lanm1_l+(k-1)] = beta;
420 T[(k-1)*lanm1_l+k] = beta;
422 int k1 = k-trlen-Ntest;
423 if ( ((k1>=0) && (k1 % cycle == 0)) || (k == lanm) || it == maxit) {
429 if (k == lanm || it == maxit) {
437 for (i=0; i<k; i++) {
438 if (Rval[i] + DBL_EPSILON >= bar) {
443 if (fabs(tr-last_tr) <= 2e-12 && jl == last_jl) {
445 fprintf(fstats,
"break: [it %d, k %d]: last_tr %.15e, tr %.15e, last_jl %d jl %d\n",
446 it, k, last_tr, tr, last_jl, jl);
451 fprintf(fstats,
"[it %d, k %d] testing: last_tr %.15e, tr %.15e, last_jl %d jl %d\n",
452 it, k, last_tr, tr, last_jl, jl);
462 for (i=0; i<k; i++) {
467 r = fabs(beta*EvecT[i*lanm1+(k-1)]);
476 fprintf(fstats,
"inner: testing conv k %4d, it %4d, count %3d jl %3d trlen %3d\n",
477 k, it, count, jl, trlen);
481 if (count*evFrac >= nev-lock) {
482 fprintf(fstats,
"count * evFrac >= nev-lock: %d * %d >= %d - %d\n", count, evFrac, nev, lock);
486 if ((count<=last_count) && (jl<=last_jl)) {
487 fprintf(fstats,
"inner: count <= last_count && jl <= last_jl: %d <= %d && %d <= %d\n", count, last_count, jl, last_jl);
505 for (i=0; i<k; i++) {
508 if (Rval[i] >= bar) {
511 DCOPY(&k, EvecT+i*lanm1_l, &one, EvecT+jl*lanm1_l, &one);
515 resi[jl] = fabs(beta*EvecT[i*lanm1_l+(k-1)]);
524 DGEMM(&cN, &cN, &n, &jl, &k, &done, V, &n, EvecT, &lanm1, &dzero, Rvec, &n);
526 DGEMM(&cN, &cN, &n, &jl, &k, &done, Z, &n, EvecT, &lanm1, &dzero, Bvec, &n);
534 for (i=0; i<jl; i++) {
536 double *q = Rvec + i*n_l;
537 double *y = Bvec + i*n_l;
544 t = sqrt(
DDOT(&n, y, &one, w2, &one));
547 t =
DNRM2(&n, y, &one);
555 DSCAL(&n, &t, y, &one);
558 DSCAL(&n, &t, w2, &one);
565 double t3 =
DDOT(&n, y, &one, w, &one);
567 if (t3 >= aa - DBL_EPSILON && t3 <= bb + DBL_EPSILON) {
573 DAXPY(&n, &nt3, w2, &one, w, &one);
576 DAXPY(&n, &nt3, y, &one, w, &one);
582 for (j=0; j<n; j++) {
588 res0 =
DNRM2(&n, w, &one);
595 nev += 1 + (int) (nev*nevInc);
597 fprintf(fstats,
"-- More eigval found: realloc space for %d evs\n", nev);
610 DCOPY(&n, y, &one, Q+lock*n_l, &one);
612 DCOPY(&n, q, &one, Y+lock*n_l, &one);
619 Rval[trlen] = Rval[i];
620 DCOPY(&n, y, &one, Z+trlen*n_l, &one);
622 DCOPY(&n, q, &one, V+trlen*n_l, &one);
626 s[trlen] = beta * EvecT[i*lanm1_l+(k-1)];
641 fprintf(fstats,
"it %4d: k %3d, jl %3d, ll %3d, lock %3d, trlen %3d\n",
642 it, k, jl, ll, lock, trlen);
645 if ((prtrlen == 0) && (ll==0)) {
648 fprintf(fstats,
"--------------------------------------------------\n");
649 fprintf(fstats,
" --> No ev.s left to be computed\n");
655 fprintf(fstats,
"--------------------------------------------------\n");
656 fprintf(fstats,
" --> Max its reached without convergence\n");
661 memset(T, 0, lanm1_l*lanm1_l*
sizeof(
double));
663 DCOPY(&n, V+k*n_l, &one, V+trlen*n_l, &one);
665 DCOPY(&n, Z+k*n_l, &one, Z+trlen*n_l, &one);
670 fprintf(fstats,
" Number of evals found = %d\n", lock);
671 fprintf(fstats,
"--------------------------------------------------\n");
void SymEigenSolver(int n, double *A, int lda, double *Q, int ldq, double *lam)
interface to LAPACK SYMMETRIC EIGEN-SOLVER
void DGEMM(char *transa, char *transb, int *m, int *n, int *k, double *alpha, double *a, int *lda, double *b, int *ldb, double *beta, double *c, int *ldc)
double DDOT(int *n, double *x, int *incx, double *y, int *incy)
void rand_double(int n, double *v)
void CGS_DGKS(int n, int k, int i_max, double *Q, double *v, double *nrmv, double *w)
Classical GS reortho with Daniel, Gragg, Kaufman, Stewart test.
void DSCAL(int *n, double *a, double *x, int *incx)
#define Calloc(base, nmem, type)
This file contains function prototypes and constant definitions internally used in EVSL...
#define Malloc(base, nmem, type)
void DGEMV(char *trans, int *m, int *n, double *alpha, double *a, int *lda, double *x, int *incx, double *beta, double *y, int *incy)
evslStat evslstat
global statistics of EVSL
void CGS_DGKS2(int n, int k, int i_max, double *Z, double *Q, double *v, double *w)
Classical GS reortho. No test. just do i_max times used in generalized ev problems.
void DAXPY(int *n, double *alpha, double *x, int *incx, double *y, int *incy)
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].
double DNRM2(int *n, double *x, int *incx)
Defs for blaslapack routines.
evslData evsldata
global variable of EVSL
#define Realloc(base, nmem, type)
void RatFiltApply(int n, ratparams *rat, double *b, double *x, double *w6)
Apply rational filter R to a vetor b.
double evsl_timer()
evsl timer for mac
void DCOPY(int *n, double *dx, int *incx, double *dy, int *incy)
parameters for rational filter