Actual source code: dsrtdf.c

slepc-3.12.1 2019-11-08
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    BDC - Block-divide and conquer (see description in README file)
 12: */

 14: #include <slepc/private/dsimpl.h>
 15: #include <slepcblaslapack.h>

 17: PetscErrorCode BDC_dsrtdf_(PetscBLASInt *k,PetscBLASInt n,PetscBLASInt n1,
 18:         PetscReal *d,PetscReal *q,PetscBLASInt ldq,PetscBLASInt *indxq,
 19:         PetscReal *rho,PetscReal *z,PetscReal *dlamda,PetscReal *w,
 20:         PetscReal *q2,PetscBLASInt *indx,PetscBLASInt *indxc,PetscBLASInt *indxp,
 21:         PetscBLASInt *coltyp,PetscReal reltol,PetscBLASInt *dz,PetscBLASInt *de,
 22:         PetscBLASInt *info)
 23: {
 24: /*  -- Routine written in LAPACK Version 3.0 style -- */
 25: /* *************************************************** */
 26: /*     Written by */
 27: /*     Michael Moldaschl and Wilfried Gansterer */
 28: /*     University of Vienna */
 29: /*     last modification: March 16, 2014 */

 31: /*     Small adaptations of original code written by */
 32: /*     Wilfried Gansterer and Bob Ward, */
 33: /*     Department of Computer Science, University of Tennessee */
 34: /*     see https://doi.org/10.1137/S1064827501399432 */
 35: /* *************************************************** */

 37: /*  Purpose */
 38: /*  ======= */

 40: /*  DSRTDF merges the two sets of eigenvalues of a rank-one modified */
 41: /*  diagonal matrix D+ZZ^T together into a single sorted set. Then it tries */
 42: /*  to deflate the size of the problem. */
 43: /*  There are two ways in which deflation can occur:  when two or more */
 44: /*  eigenvalues of D are close together or if there is a tiny entry in the */
 45: /*  Z vector.  For each such occurrence the order of the related secular */
 46: /*  equation problem is reduced by one. */

 48: /*  Arguments */
 49: /*  ========= */

 51: /*  K      (output) INTEGER */
 52: /*         The number of non-deflated eigenvalues, and the order of the */
 53: /*         related secular equation. 0 <= K <=N. */

 55: /*  N      (input) INTEGER */
 56: /*         The dimension of the diagonal matrix.  N >= 0. */

 58: /*  N1     (input) INTEGER */
 59: /*         The location of the last eigenvalue in the leading sub-matrix. */
 60: /*         min(1,N) <= N1 <= max(1,N). */

 62: /*  D      (input/output) DOUBLE PRECISION array, dimension (N) */
 63: /*         On entry, D contains the eigenvalues of the two submatrices to */
 64: /*         be combined. */
 65: /*         On exit, D contains the trailing (N-K) updated eigenvalues */
 66: /*         (those which were deflated) sorted into increasing order. */

 68: /*  Q      (input/output) DOUBLE PRECISION array, dimension (LDQ, N) */
 69: /*         On entry, Q contains the eigenvectors of two submatrices in */
 70: /*         the two square blocks with corners at (1,1), (N1,N1) */
 71: /*         and (N1+1, N1+1), (N,N). */
 72: /*         On exit, Q contains the trailing (N-K) updated eigenvectors */
 73: /*         (those which were deflated) in its last N-K columns. */

 75: /*  LDQ    (input) INTEGER */
 76: /*         The leading dimension of the array Q.  LDQ >= max(1,N). */

 78: /*  INDXQ  (input/output) INTEGER array, dimension (N) */
 79: /*         The permutation which separately sorts the two sub-problems */
 80: /*         in D into ascending order.  Note that elements in the second */
 81: /*         half of this permutation must first have N1 added to their */
 82: /*         values. Destroyed on exit. */

 84: /*  RHO    (input/output) DOUBLE PRECISION */
 85: /*         On entry, the off-diagonal element associated with the rank-1 */
 86: /*         cut which originally split the two submatrices which are now */
 87: /*         being recombined. */
 88: /*         On exit, RHO has been modified to the value required by */
 89: /*         DLAED3M (made positive and multiplied by two, such that */
 90: /*         the modification vector has norm one). */

 92: /*  Z      (input/output) DOUBLE PRECISION array, dimension (N) */
 93: /*         On entry, Z contains the updating vector. */
 94: /*         On exit, the contents of Z has been destroyed by the updating */
 95: /*         process. */

 97: /*  DLAMDA (output) DOUBLE PRECISION array, dimension (N) */
 98: /*         A copy of the first K non-deflated eigenvalues which */
 99: /*         will be used by DLAED3M to form the secular equation. */

101: /*  W      (output) DOUBLE PRECISION array, dimension (N) */
102: /*         The first K values of the final deflation-altered z-vector */
103: /*         which will be passed to DLAED3M. */

105: /*  Q2     (output) DOUBLE PRECISION array, dimension */
106: /*         ( N*N ) (if everything is deflated) or */
107: /*         ( N1*(COLTYP(1)+COLTYP(2)) + (N-N1)*(COLTYP(2)+COLTYP(3)) ) */
108: /*         (if not everything is deflated) */
109: /*         If everything is deflated, then N*N intermediate workspace */
110: /*         is needed in Q2. */
111: /*         If not everything is deflated, Q2 contains on exit a copy of the */
112: /*         first K non-deflated eigenvectors which will be used by */
113: /*         DLAED3M in a matrix multiply (DGEMM) to accumulate the new */
114: /*         eigenvectors, followed by the N-K deflated eigenvectors. */

116: /*  INDX   (workspace) INTEGER array, dimension (N) */
117: /*         The permutation used to sort the contents of DLAMDA into */
118: /*         ascending order. */

120: /*  INDXC  (output) INTEGER array, dimension (N) */
121: /*         The permutation used to arrange the columns of the deflated */
122: /*         Q matrix into three groups:  columns in the first group contain */
123: /*         non-zero elements only at and above N1 (type 1), columns in the */
124: /*         second group are dense (type 2), and columns in the third group */
125: /*         contain non-zero elements only below N1 (type 3). */

127: /*  INDXP  (workspace) INTEGER array, dimension (N) */
128: /*         The permutation used to place deflated values of D at the end */
129: /*         of the array.  INDXP(1:K) points to the nondeflated D-values */
130: /*         and INDXP(K+1:N) points to the deflated eigenvalues. */

132: /*  COLTYP (workspace/output) INTEGER array, dimension (N) */
133: /*         During execution, a label which will indicate which of the */
134: /*         following types a column in the Q2 matrix is: */
135: /*         1 : non-zero in the upper half only; */
136: /*         2 : dense; */
137: /*         3 : non-zero in the lower half only; */
138: /*         4 : deflated. */
139: /*         On exit, COLTYP(i) is the number of columns of type i, */
140: /*         for i=1 to 4 only. */

142: /*  RELTOL (input) DOUBLE PRECISION */
143: /*         User specified deflation tolerance. If RELTOL.LT.20*EPS, then */
144: /*         the standard value used in the original LAPACK routines is used. */

146: /*  DZ     (output) INTEGER, DZ.GE.0 */
147: /*         counts the deflation due to a small component in the modification */
148: /*         vector. */

150: /*  DE     (output) INTEGER, DE.GE.0 */
151: /*         counts the deflation due to close eigenvalues. */

153: /*  INFO   (output) INTEGER */
154: /*          = 0:  successful exit. */
155: /*          < 0:  if INFO = -i, the i-th argument had an illegal value. */

157: /*  Further Details */
158: /*  =============== */

160: /*  Based on code written by */
161: /*  Wilfried Gansterer and Bob Ward, */
162: /*  Department of Computer Science, University of Tennessee */

164: /*  Based on the design of the LAPACK code DLAED2 with modifications */
165: /*  to allow a block divide and conquer scheme */

167: /*  DLAED2 was written by Jeff Rutter, Computer Science Division, University */
168: /*  of California at Berkeley, USA, and modified by Francoise Tisseur, */
169: /*  University of Tennessee. */

171: /*  ===================================================================== */

173: #if defined(SLEPC_MISSING_LAPACK_LAMRG) || defined(SLEPC_MISSING_LAPACK_LACPY)
175:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAMRG/LACPY - Lapack routine is unavailable");
176: #else
177:   PetscReal    c, s, t, eps, tau, tol, dmax, dmone = -1.;
178:   PetscBLASInt i, j, i1, k2, n2, ct, nj, pj=0, js, iq1, iq2;
179:   PetscBLASInt psm[4], imax, jmax, ctot[4], factmp=1, one=1;

182:   *info = 0;

184:   if (n < 0) *info = -2;
185:   else if (n1 < PetscMin(1,n) || n1 > PetscMax(1,n)) *info = -3;
186:   else if (ldq < PetscMax(1,n)) *info = -6;
187:   if (*info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong argument %d in DSRTDF",-(*info));

189:   /* Initialize deflation counters */

191:   *dz = 0;
192:   *de = 0;

194: /* **************************************************************************** */

196:   /* Quick return if possible */

198:   if (n == 0) return(0);

200: /* **************************************************************************** */

202:   n2 = n - n1;

204:   /* Modify Z so that RHO is positive. */

206:   if (*rho < 0.) PetscStackCallBLAS("BLASscal",BLASscal_(&n2, &dmone, &z[n1], &one));

208:   /* Normalize z so that norm2(z) = 1.  Since z is the concatenation of */
209:   /* two normalized vectors, norm2(z) = sqrt(2). (NOTE that this holds also */
210:   /* here in the approximate block-tridiagonal D&C because the two vectors are */
211:   /* singular vectors, but it would not necessarily hold in a D&C for a */
212:   /* general banded matrix !) */

214:   t = 1. / PetscSqrtReal(2.);
215:   PetscStackCallBLAS("BLASscal",BLASscal_(&n, &t, z, &one));

217:   /* NOTE: at this point the value of RHO is modified in order to */
218:   /*       normalize Z:    RHO = ABS( norm2(z)**2 * RHO */
219:   /*       in our case:    norm2(z) = sqrt(2), and therefore: */

221:   *rho = PetscAbs(*rho * 2.);

223:   /* Sort the eigenvalues into increasing order */

225:   for (i = n1; i < n; ++i) indxq[i] += n1;

227:   /* re-integrate the deflated parts from the last pass */

229:   for (i = 0; i < n; ++i) dlamda[i] = d[indxq[i]-1];
230:   PetscStackCallBLAS("LAPACKlamrg",LAPACKlamrg_(&n1, &n2, dlamda, &one, &one, indxc));
231:   for (i = 0; i < n; ++i) indx[i] = indxq[indxc[i]-1];
232:   for (i = 0; i < n; ++i) indxq[i] = 0;

234:   /* Calculate the allowable deflation tolerance */

236:   /* imax = BLASamax_(&n, z, &one); */
237:   imax = 1;
238:   dmax = PetscAbsReal(z[0]);
239:   for (i=1;i<n;i++) {
240:     if (PetscAbsReal(z[i])>dmax) {
241:       imax = i+1;
242:       dmax = PetscAbsReal(z[i]);
243:     }
244:   }
245:   /* jmax = BLASamax_(&n, d, &one); */
246:   jmax = 1;
247:   dmax = PetscAbsReal(d[0]);
248:   for (i=1;i<n;i++) {
249:     if (PetscAbsReal(d[i])>dmax) {
250:       jmax = i+1;
251:       dmax = PetscAbsReal(d[i]);
252:     }
253:   }

255:   eps = LAPACKlamch_("Epsilon");
256:   if (reltol < eps * 20) {
257:     /* use the standard deflation tolerance from the original LAPACK code */
258:     tol = eps * 8. * PetscMax(PetscAbs(d[jmax-1]),PetscAbs(z[imax-1]));
259:   } else {
260:     /* otherwise set TOL to the input parameter RELTOL ! */
261:     tol = reltol * PetscMax(PetscAbs(d[jmax-1]),PetscAbs(z[imax-1]));
262:   }

264:   /* If the rank-1 modifier is small enough, no more needs to be done */
265:   /* except to reorganize Q so that its columns correspond with the */
266:   /* elements in D. */

268:   if (*rho * PetscAbs(z[imax-1]) <= tol) {

270:     /* complete deflation because of small rank-one modifier */
271:     /* NOTE: in this case we need N*N space in the array Q2 ! */

273:     *dz = n; *k = 0;
274:     iq2 = 1;
275:     for (j = 0; j < n; ++j) {
276:       i = indx[j]; indxq[j] = i;
277:       PetscStackCallBLAS("BLAScopy",BLAScopy_(&n, &q[(i-1)*ldq], &one, &q2[iq2-1], &one));
278:       dlamda[j] = d[i-1];
279:       iq2 += n;
280:     }
281:     PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("A", &n, &n, q2, &n, q, &ldq));
282:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&n, dlamda, &one, d, &one));
283:     return(0);
284:   }

286:   /* If there are multiple eigenvalues then the problem deflates.  Here */
287:   /* the number of equal eigenvalues is found.  As each equal */
288:   /* eigenvalue is found, an elementary reflector is computed to rotate */
289:   /* the corresponding eigensubspace so that the corresponding */
290:   /* components of Z are zero in this new basis. */

292:   /* initialize the column types */

294:   /* first N1 columns are initially (before deflation) of type 1 */
295:   for (i = 0; i < n1; ++i) coltyp[i] = 1;
296:   /* columns N1+1 to N are initially (before deflation) of type 3 */
297:   for (i = n1; i < n; ++i) coltyp[i] = 3;

299:   *k = 0;
300:   k2 = n + 1;
301:   for (j = 0; j < n; ++j) {
302:     nj = indx[j]-1;
303:     if (*rho * PetscAbs(z[nj]) <= tol) {

305:       /* Deflate due to small z component. */
306:       ++(*dz);
307:       --k2;

309:       /* deflated eigenpair, therefore column type 4 */
310:       coltyp[nj] = 4;
311:       indxp[k2-1] = nj+1;
312:       indxq[k2-1] = nj+1;
313:       if (j+1 == n) goto L100;
314:     } else {

316:       /* not deflated */
317:       pj = nj;
318:       goto L80;
319:     }
320:   }
321:   factmp = 1;
322: L80:
323:   ++j;
324:   nj = indx[j]-1;
325:   if (j+1 > n) goto L100;
326:   if (*rho * PetscAbs(z[nj]) <= tol) {

328:     /* Deflate due to small z component. */
329:     ++(*dz);
330:     --k2;
331:     coltyp[nj] = 4;
332:     indxp[k2-1] = nj+1;
333:     indxq[k2-1] = nj+1;
334:   } else {

336:     /* Check if eigenvalues are close enough to allow deflation. */
337:     s = z[pj];
338:     c = z[nj];

340:     /* Find sqrt(a**2+b**2) without overflow or */
341:     /* destructive underflow. */

343:     tau = LAPACKlapy2_(&c, &s);
344:     t = d[nj] - d[pj];
345:     c /= tau;
346:     s = -s / tau;
347:     if (PetscAbs(t * c * s) <= tol) {

349:       /* Deflate due to close eigenvalues. */
350:       ++(*de);
351:       z[nj] = tau;
352:       z[pj] = 0.;
353:       if (coltyp[nj] != coltyp[pj]) coltyp[nj] = 2;

355:         /* if deflation happens across the two eigenvector blocks */
356:         /* (eigenvalues corresponding to different blocks), then */
357:         /* on column in the eigenvector matrix fills up completely */
358:         /* (changes from type 1 or 3 to type 2) */

360:         /* eigenpair PJ is deflated, therefore the column type changes */
361:         /* to 4 */

363:         coltyp[pj] = 4;
364:         PetscStackCallBLAS("BLASrot",BLASrot_(&n, &q[pj*ldq], &one, &q[nj*ldq], &one, &c, &s));
365:         t = d[pj] * (c * c) + d[nj] * (s * s);
366:         d[nj] = d[pj] * (s * s) + d[nj] * (c * c);
367:         d[pj] = t;
368:         --k2;
369:         i = 1;
370: L90:
371:         if (k2 + i <= n) {
372:           if (d[pj] < d[indxp[k2 + i-1]-1]) {
373:             indxp[k2 + i - 2] = indxp[k2 + i - 1];
374:             indxq[k2 + i - 2] = indxq[k2 + i - 1];
375:             indxp[k2 + i - 1] = pj + 1;
376:             indxq[k2 + i - 2] = pj + 1;
377:             ++i;
378:             goto L90;
379:           } else {
380:             indxp[k2 + i - 2] = pj + 1;
381:             indxq[k2 + i - 2] = pj + 1;
382:           }
383:         } else {
384:           indxp[k2 + i - 2] = pj + 1;
385:           indxq[k2 + i - 2] = pj + 1;
386:         }
387:         pj = nj;
388:         factmp = -1;
389:       } else {

391:       /* do not deflate */
392:       ++(*k);
393:       dlamda[*k-1] = d[pj];
394:       w[*k-1] = z[pj];
395:       indxp[*k-1] = pj + 1;
396:       indxq[*k-1] = pj + 1;
397:       factmp = 1;
398:       pj = nj;
399:     }
400:   }
401:   goto L80;
402: L100:

404:   /* Record the last eigenvalue. */
405:   ++(*k);
406:   dlamda[*k-1] = d[pj];
407:   w[*k-1] = z[pj];
408:   indxp[*k-1] = pj+1;
409:   indxq[*k-1] = (pj+1) * factmp;

411:   /* Count up the total number of the various types of columns, then */
412:   /* form a permutation which positions the four column types into */
413:   /* four uniform groups (although one or more of these groups may be */
414:   /* empty). */

416:   for (j = 0; j < 4; ++j) ctot[j] = 0;
417:   for (j = 0; j < n; ++j) {
418:     ct = coltyp[j];
419:     ++ctot[ct-1];
420:   }

422:   /* PSM(*) = Position in SubMatrix (of types 1 through 4) */
423:   psm[0] = 1;
424:   psm[1] = ctot[0] + 1;
425:   psm[2] = psm[1] + ctot[1];
426:   psm[3] = psm[2] + ctot[2];
427:   *k = n - ctot[3];

429:   /* Fill out the INDXC array so that the permutation which it induces */
430:   /* will place all type-1 columns first, all type-2 columns next, */
431:   /* then all type-3's, and finally all type-4's. */

433:   for (j = 0; j < n; ++j) {
434:     js = indxp[j];
435:     ct = coltyp[js-1];
436:     indx[psm[ct - 1]-1] = js;
437:     indxc[psm[ct - 1]-1] = j+1;
438:     ++psm[ct - 1];
439:   }

441:   /* Sort the eigenvalues and corresponding eigenvectors into DLAMDA */
442:   /* and Q2 respectively.  The eigenvalues/vectors which were not */
443:   /* deflated go into the first K slots of DLAMDA and Q2 respectively, */
444:   /* while those which were deflated go into the last N - K slots. */

446:   i = 0;
447:   iq1 = 0;
448:   iq2 = (ctot[0] + ctot[1]) * n1;
449:   for (j = 0; j < ctot[0]; ++j) {
450:     js = indx[i];
451:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&n1, &q[(js-1)*ldq], &one, &q2[iq1], &one));
452:     z[i] = d[js-1];
453:     ++i;
454:     iq1 += n1;
455:   }

457:   for (j = 0; j < ctot[1]; ++j) {
458:     js = indx[i];
459:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&n1, &q[(js-1)*ldq], &one, &q2[iq1], &one));
460:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&n2, &q[n1+(js-1)*ldq], &one, &q2[iq2], &one));
461:     z[i] = d[js-1];
462:     ++i;
463:     iq1 += n1;
464:     iq2 += n2;
465:   }

467:   for (j = 0; j < ctot[2]; ++j) {
468:     js = indx[i];
469:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&n2, &q[n1+(js-1)*ldq], &one, &q2[iq2], &one));
470:     z[i] = d[js-1];
471:     ++i;
472:     iq2 += n2;
473:   }

475:   iq1 = iq2;
476:   for (j = 0; j < ctot[3]; ++j) {
477:     js = indx[i];
478:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&n, &q[(js-1)*ldq], &one, &q2[iq2], &one));
479:     iq2 += n;
480:     z[i] = d[js-1];
481:     ++i;
482:   }

484:   /* The deflated eigenvalues and their corresponding vectors go back */
485:   /* into the last N - K slots of D and Q respectively. */

487:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("A", &n, &ctot[3], &q2[iq1], &n, &q[*k*ldq], &ldq));
488:   i1 = n - *k;
489:   PetscStackCallBLAS("BLAScopy",BLAScopy_(&i1, &z[*k], &one, &d[*k], &one));

491:   /* Copy CTOT into COLTYP for referencing in DLAED3M. */

493:   for (j = 0; j < 4; ++j) coltyp[j] = ctot[j];
494:   return(0);
495: #endif
496: }