void rqi ( struct vtx_data **A, /* matrix/graph being analyzed */ double **yvecs, /* eigenvectors to be refined */ int index, /* index of vector in yvecs to be refined */ int n, /* number of rows/columns in matrix */ double *r1, double *r2, double *v, double *w, double *x, double *y, double *work, /* work space for symmlq */ double tol, /* error tolerance in eigenpair */ double initshift, /* initial shift */ double *evalest, /* returned eigenvalue */ double *vwsqrt, /* square roots of vertex weights */ struct orthlink *orthlist, /* lower evecs to orthogonalize against */ int cube_or_mesh, /* 0 => hypercube, d => d-dimensional mesh */ int nsets, /* number of sets to divide into */ int *assignment, /* set number of each vtx (length n+1) */ int *active, /* space for nvtxs integers */ int mediantype, /* which partitioning strategy to use */ double *goal, /* desired set sizes */ int vwgt_max, /* largest vertex weight */ int ndims /* dimensionality of partition */ ) { extern int DEBUG_EVECS; /* debug flag for eigen computation */ extern int DEBUG_TRACE; /* trace main execution path */ extern int WARNING_EVECS; /* warning flag for eigen computation */ extern int RQI_CONVERGENCE_MODE; /* type of convergence monitoring to do */ int rqisteps; /* # rqi rqisteps */ double res; /* convergence quant for rqi */ double last_res; /* res on previous rqi step */ double macheps; /* machine precision calculated by symmlq */ double normxlim; /* a stopping criteria for symmlq */ double normx; /* norm of the solution vector */ int symmlqitns; /* # symmlq itns */ int inv_it_steps; /* intial steps of inverse iteration */ long itnmin; /* symmlq input */ double shift, rtol; /* symmlq input */ long precon, goodb, nout; /* symmlq input */ long checka, intlim; /* symmlq input */ double anorm, acond; /* symmlq output */ double rnorm, ynorm; /* symmlq output */ long istop, itn; /* symmlq output */ long long_n; /* copy of n for passing to symmlq */ int warning; /* warning on possible misconvergence */ double factor; /* ratio between previous res and new tol */ double minfactor; /* minimum acceptable value of factor */ int converged; /* has process converged yet? */ double *u; /* name of vector being refined */ int *old_assignment=NULL;/* previous assignment vector */ int *assgn_pntr; /* pntr to assignment vector */ int *old_assgn_pntr; /* pntr to previous assignment vector */ int assigndiff=0; /* discrepancies between old and new assignment */ int assigntol=0; /* tolerance on convergence of assignment vector */ int first; /* is this the first RQI step? */ int i; /* loop index */ double dot(), ch_norm(); int symmlq_(); void splarax(), scadd(), vecscale(), doubleout(), assign(), x2y(), strout(); if (DEBUG_TRACE > 0) { printf("<Entering rqi>\n"); } /* Initialize RQI loop */ u = yvecs[index]; splarax(y, A, n, u, vwsqrt, r1); shift = dot(u, 1, n, y); scadd(y, 1, n, -shift, u); res = ch_norm(y, 1, n); /* eigen-residual */ rqisteps = 0; /* a counter */ symmlqitns = 0; /* a counter */ /* Set invariant symmlq parameters */ precon = FALSE; /* FALSE until we figure out a good way */ goodb = TRUE; /* should be TRUE for this application */ nout = 0; /* set to 0 for no Symmlq output; 6 for lots */ checka = FALSE; /* if don't know by now, too bad */ intlim = n; /* set to enforce a maximum number of Symmlq itns */ itnmin = 0; /* set to enforce a minimum number of Symmlq itns */ long_n = n; /* type change for alint */ if (DEBUG_EVECS > 0) { printf("Using RQI/Symmlq refinement on graph with %d vertices.\n", n); } if (DEBUG_EVECS > 1) { printf(" step lambda est. Ares Symmlq its. istop factor delta\n"); printf(" 0"); doubleout(shift, 1); doubleout(res, 1); printf("\n"); } if (RQI_CONVERGENCE_MODE == 1) { assigntol = tol * n; old_assignment = smalloc((n + 1) * sizeof(int)); } /* Perform RQI */ inv_it_steps = 2; warning = FALSE; factor = 10; minfactor = factor / 2; first = TRUE; if (res < tol) converged = TRUE; else converged = FALSE; while (!converged) { if (res / tol < 1.2) { factor = max(factor / 2, minfactor); } rtol = res / factor; /* exit Symmlq if iterate is this large */ normxlim = 1.0 / rtol; if (rqisteps < inv_it_steps) { shift = initshift; } symmlq_(&long_n, &u[1], &r1[1], &r2[1], &v[1], &w[1], &x[1], &y[1], work, &checka, &goodb, &precon, &shift, &nout, &intlim, &rtol, &istop, &itn, &anorm, &acond, &rnorm, &ynorm, (double *) A, vwsqrt, (double *) orthlist, &macheps, &normxlim, &itnmin); symmlqitns += itn; normx = ch_norm(x, 1, n); vecscale(u, 1, n, 1.0 / normx, x); splarax(y, A, n, u, vwsqrt, r1); shift = dot(u, 1, n, y); scadd(y, 1, n, -shift, u); last_res = res; res = ch_norm(y, 1, n); if (res > last_res) { warning = TRUE; } rqisteps++; if (res < tol) converged = TRUE; if (RQI_CONVERGENCE_MODE == 1 && !converged && ndims == 1) { if (first) { assign(A, yvecs, n, 1, cube_or_mesh, nsets, vwsqrt, assignment, active, mediantype, goal, vwgt_max); x2y(yvecs, ndims, n, vwsqrt); first = FALSE; assigndiff = n; /* dummy value for debug chart */ } else { /* copy assignment to old_assignment */ assgn_pntr = assignment; old_assgn_pntr = old_assignment; for (i = n + 1; i; i--) { *old_assgn_pntr++ = *assgn_pntr++; } assign(A, yvecs, n, ndims, cube_or_mesh, nsets, vwsqrt, assignment, active, mediantype, goal, vwgt_max); x2y(yvecs, ndims, n, vwsqrt); /* count differences in assignment */ assigndiff = 0; assgn_pntr = assignment; old_assgn_pntr = old_assignment; for (i = n + 1; i; i--) { if (*old_assgn_pntr++ != *assgn_pntr++) assigndiff++; } assigndiff = min(assigndiff, n - assigndiff); if (assigndiff <= assigntol) converged = TRUE; } } if (DEBUG_EVECS > 1) { printf(" %2d", rqisteps); doubleout(shift, 1); doubleout(res, 1); printf(" %3ld", itn); printf(" %ld", istop); printf(" %g", factor); if (RQI_CONVERGENCE_MODE == 1) printf(" %d\n", assigndiff); else printf("\n"); } } *evalest = shift; if (WARNING_EVECS > 0 && warning) { strout("WARNING: Residual convergence not monotonic; RQI may have misconverged.\n"); } if (DEBUG_EVECS > 0) { printf("Eval "); doubleout(*evalest, 1); printf(" RQI steps %d, Symmlq iterations %d.\n\n", rqisteps, symmlqitns); } if (RQI_CONVERGENCE_MODE == 1) { sfree(old_assignment); } }
void lanczos_FO ( struct vtx_data **A, /* graph data structure */ int n, /* number of rows/colums in matrix */ int d, /* problem dimension = # evecs to find */ double **y, /* columns of y are eigenvectors of A */ double *lambda, /* ritz approximation to eigenvals of A */ double *bound, /* on ritz pair approximations to eig pairs of A */ double eigtol, /* tolerance on eigenvectors */ double *vwsqrt, /* square root of vertex weights */ double maxdeg, /* maximum degree of graph */ int version /* 1 = standard mode, 2 = inverse operator mode */ ) { extern FILE *Output_File; /* output file or NULL */ extern int DEBUG_EVECS; /* print debugging output? */ extern int DEBUG_TRACE; /* trace main execution path */ extern int WARNING_EVECS; /* print warning messages? */ extern int LANCZOS_MAXITNS; /* maximum Lanczos iterations allowed */ extern double BISECTION_SAFETY; /* safety factor for bisection algorithm */ extern double SRESTOL; /* resid tol for T evec comp */ extern double DOUBLE_MAX; /* Warning on inaccurate computation of evec of T */ extern double splarax_time; /* time matvecs */ extern double orthog_time; /* time orthogonalization work */ extern double tevec_time; /* time tridiagonal eigvec work */ extern double evec_time; /* time to generate eigenvectors */ extern double ql_time; /* time tridiagonal eigval work */ extern double blas_time; /* time for blas (not assembly coded) */ extern double init_time; /* time for allocating memory, etc. */ extern double scan_time; /* time for scanning bounds list */ extern double debug_time; /* time for debug computations and output */ int i, j; /* indicies */ int maxj; /* maximum number of Lanczos iterations */ double *u, *r; /* Lanczos vectors */ double *Aq; /* sparse matrix-vector product vector */ double *alpha, *beta; /* the Lanczos scalars from each step */ double *ritz; /* copy of alpha for tqli */ double *workj; /* work vector (eg. for tqli) */ double *workn; /* work vector (eg. for checkeig) */ double *s; /* eigenvector of T */ double **q; /* columns of q = Lanczos basis vectors */ double *bj; /* beta(j)*(last element of evecs of T) */ double bis_safety; /* real safety factor for bisection algorithm */ double Sres; /* how well Tevec calculated eigvecs */ double Sres_max; /* Maximum value of Sres */ int inc_bis_safety; /* need to increase bisection safety */ double *Ares; /* how well Lanczos calculated each eigpair */ double *inv_lambda; /* eigenvalues of inverse operator */ int *index; /* the Ritz index of an eigenpair */ struct orthlink *orthlist = NULL; /* vectors to orthogonalize against in Lanczos */ struct orthlink *orthlist2 = NULL; /* vectors to orthogonalize against in Symmlq */ struct orthlink *temp; /* for expanding orthogonalization list */ double *ritzvec=NULL; /* ritz vector for current iteration */ double *zeros=NULL; /* vector of all zeros */ double *ones=NULL; /* vector of all ones */ struct scanlink *scanlist; /* list of fields for min ritz vals */ struct scanlink *curlnk; /* for traversing the scanlist */ double bji_tol; /* tol on bji estimate of A e-residual */ int converged; /* has the iteration converged? */ double time; /* current clock time */ double shift, rtol; /* symmlq input */ long precon, goodb, nout; /* symmlq input */ long checka, intlim; /* symmlq input */ double anorm, acond; /* symmlq output */ double rnorm, ynorm; /* symmlq output */ long istop, itn; /* symmlq output */ double macheps; /* machine precision calculated by symmlq */ double normxlim; /* a stopping criteria for symmlq */ long itnmin; /* enforce minimum number of iterations */ int symmlqitns; /* # symmlq itns */ double *wv1=NULL, *wv2=NULL, *wv3=NULL; /* Symmlq work space */ double *wv4=NULL, *wv5=NULL, *wv6=NULL; /* Symmlq work space */ long long_n; /* long int copy of n for symmlq */ int ritzval_flag = 0; /* status flag for ql() */ double Anorm; /* Norm estimate of the Laplacian matrix */ int left, right; /* ranges on the search for ritzvals */ int memory_ok; /* TRUE as long as don't run out of memory */ double *mkvec(); /* allocates space for a vector */ double *mkvec_ret(); /* mkvec() which returns error code */ double dot(); /* standard dot product routine */ struct orthlink *makeorthlnk(); /* make space for entry in orthog. set */ double ch_norm(); /* vector norm */ double Tevec(); /* calc evec of T by linear recurrence */ struct scanlink *mkscanlist(); /* make scan list for min ritz vecs */ double lanc_seconds(); /* current clock timer */ int symmlq_(), get_ritzvals(); void setvec(), vecscale(), update(), vecran(), strout(); void splarax(), scanmin(), scanmax(), frvec(), orthogonalize(); void orthog1(), orthogvec(), bail(), warnings(), mkeigvecs(); if (DEBUG_TRACE > 0) { printf("<Entering lanczos_FO>\n"); } if (DEBUG_EVECS > 0) { if (version == 1) { printf("Full orthogonalization Lanczos, matrix size = %d\n", n); } else { printf("Full orthogonalization Lanczos, inverted operator, matrix size = %d\n", n); } } /* Initialize time. */ time = lanc_seconds(); if (n < d + 1) { bail("ERROR: System too small for number of eigenvalues requested.",1); /* d+1 since don't use zero eigenvalue pair */ } /* Allocate Lanczos space. */ maxj = LANCZOS_MAXITNS; u = mkvec(1, n); r = mkvec(1, n); Aq = mkvec(1, n); ritzvec = mkvec(1, n); zeros = mkvec(1, n); setvec(zeros, 1, n, 0.0); workn = mkvec(1, n); Ares = mkvec(1, d); inv_lambda = mkvec(1, d); index = smalloc((d + 1) * sizeof(int)); alpha = mkvec(1, maxj); beta = mkvec(1, maxj + 1); ritz = mkvec(1, maxj); s = mkvec(1, maxj); bj = mkvec(1, maxj); workj = mkvec(1, maxj + 1); q = smalloc((maxj + 1) * sizeof(double *)); scanlist = mkscanlist(d); if (version == 2) { /* Allocate Symmlq space all in one chunk. */ wv1 = smalloc(6 * (n + 1) * sizeof(double)); wv2 = &wv1[(n + 1)]; wv3 = &wv1[2 * (n + 1)]; wv4 = &wv1[3 * (n + 1)]; wv5 = &wv1[4 * (n + 1)]; wv6 = &wv1[5 * (n + 1)]; /* Set invariant symmlq parameters */ precon = FALSE; /* FALSE until we figure out a good way */ goodb = FALSE; /* should be FALSE for this application */ checka = FALSE; /* if don't know by now, too bad */ intlim = n; /* set to enforce a maximum number of Symmlq itns */ itnmin = 0; /* set to enforce a minimum number of Symmlq itns */ shift = 0.0; /* since just solving rather than doing RQI */ symmlqitns = 0; /* total number of Symmlq iterations */ nout = 0; /* Effectively disabled - see notes in symmlq.f */ rtol = 1.0e-5; /* requested residual tolerance */ normxlim = DOUBLE_MAX; /* Effectively disables ||x|| termination criterion */ long_n = n; /* copy to long for linting */ } /* Initialize. */ vecran(r, 1, n); if (vwsqrt == NULL) { /* whack one's direction from initial vector */ orthog1(r, 1, n); /* list the ones direction for later use in Symmlq */ if (version == 2) { orthlist2 = makeorthlnk(); ones = mkvec(1, n); setvec(ones, 1, n, 1.0); orthlist2->vec = ones; orthlist2->pntr = NULL; } } else { /* whack vwsqrt direction from initial vector */ orthogvec(r, 1, n, vwsqrt); if (version == 2) { /* list the vwsqrt direction for later use in Symmlq */ orthlist2 = makeorthlnk(); orthlist2->vec = vwsqrt; orthlist2->pntr = NULL; } } beta[1] = ch_norm(r, 1, n); q[0] = zeros; bji_tol = eigtol; orthlist = NULL; Sres_max = 0.0; Anorm = 2 * maxdeg; /* Gershgorin estimate for ||A|| */ bis_safety = BISECTION_SAFETY; inc_bis_safety = FALSE; init_time += lanc_seconds() - time; /* Main Lanczos loop. */ j = 1; converged = FALSE; memory_ok = TRUE; while ((j <= maxj) && (converged == FALSE) && memory_ok) { time = lanc_seconds(); /* Allocate next Lanczos vector. If fail, back up one step and compute approx. eigvec. */ q[j] = mkvec_ret(1, n); if (q[j] == NULL) { memory_ok = FALSE; if (DEBUG_EVECS > 0 || WARNING_EVECS > 0) { strout("WARNING: Lanczos out of memory; computing best approximation available.\n"); } if (j <= 2) { bail("ERROR: Sorry, can't salvage Lanczos.",1); /* ... save yourselves, men. */ } j--; } vecscale(q[j], 1, n, 1.0 / beta[j], r); blas_time += lanc_seconds() - time; time = lanc_seconds(); if (version == 1) { splarax(Aq, A, n, q[j], vwsqrt, workn); } else { symmlq_(&long_n, &(q[j][1]), &wv1[1], &wv2[1], &wv3[1], &wv4[1], &Aq[1], &wv5[1], &wv6[1], &checka, &goodb, &precon, &shift, &nout, &intlim, &rtol, &istop, &itn, &anorm, &acond, &rnorm, &ynorm, (double *) A, vwsqrt, (double *) orthlist2, &macheps, &normxlim, &itnmin); symmlqitns += itn; if (DEBUG_EVECS > 2) { printf("Symmlq report: rtol %g\n", rtol); printf(" system norm %g, solution norm %g\n", anorm, ynorm); printf(" system condition %g, residual %g\n", acond, rnorm); printf(" termination condition %2ld, iterations %3ld\n", istop, itn); } } splarax_time += lanc_seconds() - time; time = lanc_seconds(); update(u, 1, n, Aq, -beta[j], q[j - 1]); alpha[j] = dot(u, 1, n, q[j]); update(r, 1, n, u, -alpha[j], q[j]); blas_time += lanc_seconds() - time; time = lanc_seconds(); if (vwsqrt == NULL) { orthog1(r, 1, n); } else { orthogvec(r, 1, n, vwsqrt); } orthogonalize(r, n, orthlist); temp = orthlist; orthlist = makeorthlnk(); orthlist->vec = q[j]; orthlist->pntr = temp; beta[j + 1] = ch_norm(r, 1, n); orthog_time += lanc_seconds() - time; time = lanc_seconds(); left = j/2; right = j - left + 1; if (inc_bis_safety) { bis_safety *= 10; inc_bis_safety = FALSE; } ritzval_flag = get_ritzvals(alpha, beta+1, j, Anorm, workj+1, ritz, d, left, right, eigtol, bis_safety); /* ... have to off-set beta and workj since full orthogonalization indexes these from 1 to maxj+1 whereas selective orthog. indexes them from 0 to maxj */ if (ritzval_flag != 0) { bail("ERROR: Both Sturm bisection and QL failed.",1); /* ... give up. */ } ql_time += lanc_seconds() - time; /* Convergence check using Paige bji estimates. */ time = lanc_seconds(); for (i = 1; i <= j; i++) { Sres = Tevec(alpha, beta, j, ritz[i], s); if (Sres > Sres_max) { Sres_max = Sres; } if (Sres > SRESTOL) { inc_bis_safety = TRUE; } bj[i] = s[j] * beta[j + 1]; } tevec_time += lanc_seconds() - time; time = lanc_seconds(); if (version == 1) { scanmin(ritz, 1, j, &scanlist); } else { scanmax(ritz, 1, j, &scanlist); } converged = TRUE; if (j < d) converged = FALSE; else { curlnk = scanlist; while (curlnk != NULL) { if (bj[curlnk->indx] > bji_tol) { converged = FALSE; } curlnk = curlnk->pntr; } } scan_time += lanc_seconds() - time; j++; } j--; /* Collect eigenvalue and bound information. */ time = lanc_seconds(); mkeigvecs(scanlist,lambda,bound,index,bj,d,&Sres_max,alpha,beta+1,j,s,y,n,q); evec_time += lanc_seconds() - time; /* Analyze computation for and report additional problems */ time = lanc_seconds(); if (DEBUG_EVECS>0 && version == 2) { printf("\nTotal Symmlq iterations %3d\n", symmlqitns); } if (version == 2) { for (i = 1; i <= d; i++) { lambda[i] = 1.0/lambda[i]; } } warnings(workn, A, y, n, lambda, vwsqrt, Ares, bound, index, d, j, maxj, Sres_max, eigtol, u, Anorm, Output_File); debug_time += lanc_seconds() - time; /* Free any memory allocated in this routine. */ time = lanc_seconds(); frvec(u, 1); frvec(r, 1); frvec(Aq, 1); frvec(ritzvec, 1); frvec(zeros, 1); if (vwsqrt == NULL && version == 2) { frvec(ones, 1); } frvec(workn, 1); frvec(Ares, 1); frvec(inv_lambda, 1); sfree(index); frvec(alpha, 1); frvec(beta, 1); frvec(ritz, 1); frvec(s, 1); frvec(bj, 1); frvec(workj, 1); if (version == 2) { frvec(wv1, 0); } while (scanlist != NULL) { curlnk = scanlist->pntr; sfree(scanlist); scanlist = curlnk; } for (i = 1; i <= j; i++) { frvec(q[i], 1); } while (orthlist != NULL) { temp = orthlist->pntr; sfree(orthlist); orthlist = temp; } while (version == 2 && orthlist2 != NULL) { temp = orthlist2->pntr; sfree(orthlist2); orthlist2 = temp; } sfree(q); init_time += lanc_seconds() - time; }