Actual source code: dibtdc.c

slepc-3.11.2 2019-07-30
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: static PetscErrorCode cutlr_(PetscBLASInt start,PetscBLASInt n,PetscBLASInt blkct,
 18:         PetscBLASInt *bsizes,PetscBLASInt *ranks,PetscBLASInt *cut,
 19:         PetscBLASInt *lsum,PetscBLASInt *lblks,PetscBLASInt *info)
 20: {
 21: /*  -- Routine written in LAPACK Version 3.0 style -- */
 22: /* *************************************************** */
 23: /*     Written by */
 24: /*     Michael Moldaschl and Wilfried Gansterer */
 25: /*     University of Vienna */
 26: /*     last modification: March 16, 2014 */

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

 34: /*  Purpose */
 35: /*  ======= */

 37: /*  CUTLR computes the optimal cut in a sequence of BLKCT neighboring */
 38: /*  blocks whose sizes are given by the array BSIZES. */
 39: /*  The sum of all block sizes in the sequence considered is given by N. */
 40: /*  The cut is optimal in the sense that the difference of the sizes of */
 41: /*  the resulting two halves is minimum over all cuts with minimum ranks */
 42: /*  between blocks of the sequence considered. */

 44: /*  Arguments */
 45: /*  ========= */

 47: /*  START  (input) INTEGER */
 48: /*         In the original array KSIZES of the calling routine DIBTDC, */
 49: /*         the position where the sequence considered in this routine starts. */
 50: /*         START >= 1. */

 52: /*  N      (input) INTEGER */
 53: /*         The sum of all the block sizes of the sequence to be cut = */
 54: /*         = sum_{i=1}^{BLKCT} BSIZES( I ). */
 55: /*         N >= 3. */

 57: /*  BLKCT  (input) INTEGER */
 58: /*         The number of blocks in the sequence to be cut. */
 59: /*         BLKCT >= 3. */

 61: /*  BSIZES (input) INTEGER array, dimension (BLKCT) */
 62: /*         The dimensions of the (quadratic) blocks of the sequence to be */
 63: /*         cut. sum_{i=1}^{BLKCT} BSIZES( I ) = N. */

 65: /*  RANKS  (input) INTEGER array, dimension (BLKCT-1) */
 66: /*         The ranks determining the approximations of the off-diagonal */
 67: /*         blocks in the sequence considered. */

 69: /*  CUT    (output) INTEGER */
 70: /*         After the optimum cut has been determined, the position (in the */
 71: /*         overall problem as worked on in DIBTDC !) of the last block in */
 72: /*         the first half of the sequence to be cut. */
 73: /*         START <= CUT <= START+BLKCT-2. */

 75: /*  LSUM   (output) INTEGER */
 76: /*         After the optimum cut has been determined, the sum of the */
 77: /*         block sizes in the first half of the sequence to be cut. */
 78: /*         LSUM < N. */

 80: /*  LBLKS  (output) INTEGER */
 81: /*         After the optimum cut has been determined, the number of the */
 82: /*         blocks in the first half of the sequence to be cut. */
 83: /*         1 <= LBLKS < BLKCT. */

 85: /*  INFO   (output) INTEGER */
 86: /*          = 0:  successful exit. */
 87: /*          < 0:  illegal arguments. */
 88: /*                if INFO = -i, the i-th (input) argument had an illegal */
 89: /*                value. */
 90: /*          > 0:  illegal results. */
 91: /*                if INFO = i, the i-th (output) argument had an illegal */
 92: /*                value. */

 94: /*  Further Details */
 95: /*  =============== */

 97: /*  Based on code written by */
 98: /*     Wilfried Gansterer and Bob Ward, */
 99: /*     Department of Computer Science, University of Tennessee */

101: /*  ===================================================================== */

103:   PetscBLASInt i, ksk, kchk, ksum, nhalf, deviat, mindev, minrnk, tmpsum;

106:   *info = 0;
107:   *lblks = 1;
108:   *lsum = 1;
109:   *cut = start;

111:   if (start < 1) {
112:     *info = -1;
113:   } else if (n < 3) {
114:     *info = -2;
115:   } else if (blkct < 3) {
116:     *info = -3;
117:   }
118:   if (*info == 0) {
119:     ksum = 0;
120:     kchk = 0;
121:     for (i = 0; i < blkct; ++i) {
122:       ksk = bsizes[i];
123:       ksum += ksk;
124:       if (ksk < 1) kchk = 1;
125:     }
126:     if (ksum != n || kchk == 1) *info = -4;
127:   }
128:   if (*info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong argument %d in CUTLR",-(*info));

130:   /* determine smallest rank in the range considered */

132:   minrnk = n;
133:   for (i = 0; i < blkct-1; ++i) {
134:     if (ranks[i] < minrnk) minrnk = ranks[i];
135:   }

137:   /* determine best cut among those with smallest rank */

139:   nhalf = n / 2;
140:   tmpsum = 0;
141:   mindev = n;
142:   for (i = 0; i < blkct; ++i) {
143:     tmpsum += bsizes[i];
144:     if (ranks[i] == minrnk) {

146:       /* determine deviation from "optimal" cut NHALF */

148:       deviat = tmpsum - nhalf;
149:       if (deviat<0) deviat = -deviat;

151:       /* compare to best deviation so far */

153:       if (deviat < mindev) {
154:         mindev = deviat;
155:         *cut = start + i;
156:         *lblks = i + 1;
157:         *lsum = tmpsum;
158:       }
159:     }
160:   }

162:   if (*cut < start || *cut >= start + blkct - 1) {
163:     *info = 6;
164:   } else if (*lsum < 1 || *lsum >= n) {
165:     *info = 7;
166:   } else if (*lblks < 1 || *lblks >= blkct) {
167:     *info = 8;
168:   }
169:   return(0);
170: }

172: PetscErrorCode BDC_dibtdc_(const char *jobz,PetscBLASInt n,PetscBLASInt nblks,
173:         PetscBLASInt *ksizes,PetscReal *d,PetscBLASInt l1d,PetscBLASInt l2d,
174:         PetscReal *e,PetscBLASInt *rank,PetscBLASInt l1e,PetscBLASInt l2e,
175:         PetscReal tol,PetscReal *ev,PetscReal *z,PetscBLASInt ldz,PetscReal *work,
176:         PetscBLASInt lwork,PetscBLASInt *iwork,PetscBLASInt liwork,
177:         PetscBLASInt *info,PetscBLASInt jobz_len)
178: {
179: /*  -- Routine written in LAPACK Version 3.0 style -- */
180: /* *************************************************** */
181: /*     Written by */
182: /*     Michael Moldaschl and Wilfried Gansterer */
183: /*     University of Vienna */
184: /*     last modification: March 16, 2014 */

186: /*     Small adaptations of original code written by */
187: /*     Wilfried Gansterer and Bob Ward, */
188: /*     Department of Computer Science, University of Tennessee */
189: /*     see https://doi.org/10.1137/S1064827501399432 */
190: /* *************************************************** */

192: /*  Purpose */
193: /*  ======= */

195: /*  DIBTDC computes all eigenvalues and corresponding eigenvectors of a */
196: /*  symmetric irreducible block tridiagonal matrix with rank RANK matrices */
197: /*  as the subdiagonal blocks using a block divide and conquer method. */

199: /*  Arguments */
200: /*  ========= */

202: /*  JOBZ    (input) CHARACTER*1 */
203: /*          = 'N':  Compute eigenvalues only (not implemented); */
204: /*          = 'D':  Compute eigenvalues and eigenvectors. */
205: /*                  Eigenvectors are accumulated in the */
206: /*                  divide-and-conquer process. */

208: /*  N      (input) INTEGER */
209: /*         The dimension of the symmetric irreducible block tridiagonal */
210: /*         matrix.  N >= 2. */

212: /*  NBLKS  (input) INTEGER, 2 <= NBLKS <= N */
213: /*         The number of diagonal blocks in the matrix. */

215: /*  KSIZES (input) INTEGER array, dimension (NBLKS) */
216: /*         The dimension of the square diagonal blocks from top left */
217: /*         to bottom right.  KSIZES(I) >= 1 for all I, and the sum of */
218: /*         KSIZES(I) for I = 1 to NBLKS has to be equal to N. */

220: /*  D      (input) DOUBLE PRECISION array, dimension (L1D,L2D,NBLKS) */
221: /*         The lower triangular elements of the symmetric diagonal */
222: /*         blocks of the block tridiagonal matrix.  Elements of the top */
223: /*         left diagonal block, which is of dimension KSIZES(1), are */
224: /*         contained in D(*,*,1); the elements of the next diagonal */
225: /*         block, which is of dimension KSIZES(2), are contained in */
226: /*         D(*,*,2); etc. */

228: /*  L1D    (input) INTEGER */
229: /*         The leading dimension of the array D.  L1D >= max(3,KMAX), */
230: /*         where KMAX is the dimension of the largest diagonal block. */

232: /*  L2D    (input) INTEGER */
233: /*         The second dimension of the array D.  L2D >= max(3,KMAX), */
234: /*         where KMAX is as stated in L1D above. */

236: /*  E      (input) DOUBLE PRECISION array, dimension (L1E,L2E,NBLKS-1) */
237: /*         Contains the elements of the scalars (singular values) and */
238: /*         vectors (singular vectors) defining the rank RANK subdiagonal */
239: /*         blocks of the matrix. */
240: /*         E(1:RANK(K),RANK(K)+1,K) holds the RANK(K) scalars, */
241: /*         E(:,1:RANK(K),K) holds the RANK(K) column vectors, and */
242: /*         E(:,RANK(K)+2:2*RANK(K)+1,K) holds the row vectors for the K-th */
243: /*         subdiagonal block. */

245: /*  RANK   (input) INTEGER array, dimension (NBLKS-1). */
246: /*         The ranks of all the subdiagonal blocks contained in the array E. */
247: /*         RANK( K ) <= MIN( KSIZES( K ), KSIZES( K+1 ) ) */

249: /*  L1E    (input) INTEGER */
250: /*         The leading dimension of the array E.  L1E >= max(3,2*KMAX+1), */
251: /*         where KMAX is as stated in L1D above. */

253: /*  L2E    (input) INTEGER */
254: /*         The second dimension of the array E.  L2E >= max(3,2*KMAX+1), */
255: /*         where KMAX is as stated in L1D above. */

257: /*  TOL    (input) DOUBLE PRECISION, TOL <= 1.0D-1 */
258: /*         User specified deflation tolerance for the routine DMERG2. */
259: /*         If ( 1.0D-1 >= TOL >= 20*EPS ) then TOL is used as */
260: /*         the deflation tolerance in DSRTDF. */
261: /*         If ( TOL < 20*EPS ) then the standard deflation tolerance from */
262: /*         LAPACK is used as the deflation tolerance in DSRTDF. */

264: /*  EV     (output) DOUBLE PRECISION array, dimension (N) */
265: /*         If INFO = 0, then EV contains the eigenvalues of the */
266: /*         symmetric block tridiagonal matrix in ascending order. */

268: /*  Z      (input/output) DOUBLE PRECISION array, dimension (LDZ, N) */
269: /*         On entry, Z will be the identity matrix. */
270: /*         On exit, Z contains the eigenvectors of the block tridiagonal */
271: /*         matrix. */

273: /*  LDZ    (input) INTEGER */
274: /*         The leading dimension of the array Z.  LDZ >= max(1,N). */

276: /*  WORK   (workspace) DOUBLE PRECISION array, dimension (LWORK) */

278: /*  LWORK   (input) INTEGER */
279: /*          The dimension of the array WORK. */
280: /*          In order to guarantee correct results in all cases, */
281: /*          LWORK must be at least ( 2*N**2 + 3*N ). In many cases, */
282: /*          less workspace is required. The absolute minimum required is */
283: /*          ( N**2 + 3*N ). */
284: /*          If the workspace provided is not sufficient, the routine will */
285: /*          return a corresponding error code and report how much workspace */
286: /*          was missing (see INFO). */

288: /*  IWORK  (workspace) INTEGER array, dimension (LIWORK) */

290: /*  LIWORK  (input) INTEGER */
291: /*          The dimension of the array IWORK. */
292: /*          LIWORK must be at least ( 5*N + 3 + 4*NBLKS - 4 ): */
293: /*                 5*KMAX+3 for DSYEVD, 5*N for ????, */
294: /*                 4*NBLKS-4 for the preprocessing (merging order) */
295: /*          Summarizing, the minimum integer workspace needed is */
296: /*          MAX( 5*N, 5*KMAX + 3 ) + 4*NBLKS - 4 */

298: /*  INFO   (output) INTEGER */
299: /*          = 0:  successful exit. */
300: /*          < 0, > -99: illegal arguments. */
301: /*                if INFO = -i, the i-th argument had an illegal value. */
302: /*          = -99: error in the preprocessing (call of CUTLR). */
303: /*          < -200: not enough workspace. Space for ABS(INFO + 200) */
304: /*                numbers is required in addition to the workspace provided, */
305: /*                otherwise some eigenvectors will be incorrect. */
306: /*          > 0:  The algorithm failed to compute an eigenvalue while */
307: /*                working on the submatrix lying in rows and columns */
308: /*                INFO/(N+1) through mod(INFO,N+1). */

310: /*  Further Details */
311: /*  =============== */

313: /*  Based on code written by */
314: /*     Wilfried Gansterer and Bob Ward, */
315: /*     Department of Computer Science, University of Tennessee */

317: /*  This routine is comparable to Dlaed0.f from LAPACK. */

319: /*  ===================================================================== */

321: #if defined(SLEPC_MISSING_LAPACK_LACPY) || defined(PETSC_MISSING_LAPACK_SYEV)
323:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LACPY/SYEV - Lapack routine is unavailable");
324: #else
325:   PetscBLASInt   i, j, k, np, rp1, ksk, one=1;
326:   PetscBLASInt   cut, mat1, kchk, kbrk, blks, kmax, icut, size, ksum, lsum;
327:   PetscBLASInt   lblks, rblks, isize, lwmin, ilsum;
328:   PetscBLASInt   start, vstrt, istck1, istck2, istck3, merged;
329:   PetscBLASInt   liwmin, matsiz, startp, istrtp;
330:   PetscReal      rho, done=1.0, dmone=-1.0;

334:   *info = 0;

336:   if (*(unsigned char *)jobz != 'N' && *(unsigned char *)jobz != 'D') {
337:     *info = -1;
338:   } else if (n < 2) {
339:     *info = -2;
340:   } else if (nblks < 2 || nblks > n) {
341:     *info = -3;
342:   }
343:   if (*info == 0) {
344:     ksum = 0;
345:     kmax = 0;
346:     kchk = 0;
347:     for (k = 0; k < nblks; ++k) {
348:       ksk = ksizes[k];
349:       ksum += ksk;
350:       if (ksk > kmax) kmax = ksk;
351:       if (ksk < 1) kchk = 1;
352:     }
353:     lwmin = n*n + n * 3;
354:     liwmin = PetscMax(n * 5,kmax * 5 + 3) + 4*nblks - 4;
355:     if (ksum != n || kchk == 1) {
356:       *info = -4;
357:     } else if (l1d < PetscMax(3,kmax)) {
358:       *info = -6;
359:     } else if (l2d < PetscMax(3,kmax)) {
360:       *info = -7;
361:     } else if (l1e < PetscMax(3,2*kmax + 1)) {
362:       *info = -10;
363:     } else if (l2e < PetscMax(3,2*kmax + 1)) {
364:       *info = -11;
365:     } else if (tol > .1) {
366:       *info = -12;
367:     } else if (ldz < PetscMax(1,n)) {
368:       *info = -15;
369:     } else if (lwork < lwmin) {
370:       *info = -17;
371:     } else if (liwork < liwmin) {
372:       *info = -19;
373:     }
374:   }
375:   if (*info == 0) {
376:     for (k = 0; k < nblks-1; ++k) {
377:       if (rank[k] > PetscMin(ksizes[k],ksizes[k+1]) || rank[k] < 1) *info = -9;
378:     }
379:   }

381:   if (*info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong argument %d in DIBTDC",-(*info));

383: /* **************************************************************************** */

385:   /* ...Preprocessing..................................................... */
386:   /*    Determine the optimal order for merging the subblocks and how much */
387:   /*    workspace will be needed for the merging (determined by the last */
388:   /*    merge). Cutpoints for the merging operations are determined and stored */
389:   /*    in reverse chronological order (starting with the final merging */
390:   /*    operation). */

392:   /*    integer workspace requirements for the preprocessing: */
393:   /*         4*(NBLKS-1) for merging history */
394:   /*         at most 3*(NBLKS-1) for stack */

396:   start = 1;
397:   size = n;
398:   blks = nblks;
399:   merged = 0;
400:   k = 0;

402:   /* integer workspace used for the stack is not needed any more after the */
403:   /* preprocessing and therefore can use part of the 5*N */
404:   /* integer workspace needed later on in the code */

406:   istck1 = 0;
407:   istck2 = istck1 + nblks;
408:   istck3 = istck2 + nblks;

410:   /* integer workspace used for storing the order of merges starts AFTER */
411:   /* the integer workspace 5*N+3 which is needed later on in the code */
412:   /* (5*KMAX+3 for DSYEVD, 4*N in DMERG2) */

414:   istrtp = n * 5 + 4;
415:   icut = istrtp + nblks - 1;
416:   isize = icut + nblks - 1;
417:   ilsum = isize + nblks - 1;

419: L200:

421:   if (nblks >= 3) {

423:     /* Determine the cut point. Note that in the routine CUTLR it is */
424:     /* chosen such that it yields the best balanced merging operation */
425:     /* among all the rank modifications with minimum rank. */

427:     cutlr_(start, size, blks, &ksizes[start-1], &rank[start-1], &cut,
428:                   &lsum, &lblks, info);
429:     if (*info) SETERRQ1(PETSC_COMM_SELF,1,"dibtdc: Error in cutlr, info = %d",*info);

431:   } else {
432:     cut = 1;
433:     lsum = ksizes[0];
434:     lblks = 1;
435:   }

437:   ++merged;
438:   startp = 0;
439:   for (i = 0; i < start-1; ++i) startp += ksizes[i];
440:   iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
441:   iwork[icut + (nblks - 1) - merged-1] = cut;
442:   iwork[isize + (nblks - 1) - merged-1] = size;
443:   iwork[ilsum + (nblks - 1) - merged-1] = lsum;

445:   if (lblks == 2) {

447:     /* one merge in left branch, left branch done */
448:     ++merged;
449:     iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
450:     iwork[icut + (nblks - 1) - merged-1] = start;
451:     iwork[isize + (nblks - 1) - merged-1] = lsum;
452:     iwork[ilsum + (nblks - 1) - merged-1] = ksizes[start-1];
453:   }

455:   if (lblks == 1 || lblks == 2) {

457:     /* left branch done, continue on the right side */
458:     start += lblks;
459:     size -= lsum;
460:     blks -= lblks;

462:     if (blks <= 0) SETERRQ1(PETSC_COMM_SELF,1,"dibtdc: Error in preprocessing, blks = %d",blks);

464:     if (blks == 2) {

466:       /* one merge in right branch, right branch done */
467:       ++merged;
468:       startp += lsum;
469:       iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
470:       iwork[icut + (nblks - 1) - merged-1] = start;
471:       iwork[isize + (nblks - 1) - merged-1] = size;
472:       iwork[ilsum + (nblks - 1) - merged-1] = ksizes[start-1];
473:     }

475:     if (blks == 1 || blks == 2) {

477:       /* get the next subproblem from the stack or finished */

479:       if (k >= 1) {

481:         /* something left on the stack */
482:         start = iwork[istck1 + k-1];
483:         size = iwork[istck2 + k-1];
484:         blks = iwork[istck3 + k-1];
485:         --k;
486:         goto L200;
487:       } else {

489:         /* nothing left on the stack */
490:         if (merged != nblks-1) SETERRQ(PETSC_COMM_SELF,1,"ERROR in preprocessing - not enough merges performed");

492:         /* exit preprocessing */

494:       }
495:     } else {

497:       /* BLKS.GE.3, and therefore analyze the right side */

499:       goto L200;
500:     }
501:   } else {

503:     /* LBLKS.GE.3, and therefore check the right side and */
504:     /* put it on the stack if required */

506:     rblks = blks - lblks;
507:     if (rblks >= 3) {
508:       ++k;
509:       iwork[istck1 + k-1] = cut + 1;
510:       iwork[istck2 + k-1] = size - lsum;
511:       iwork[istck3 + k-1] = rblks;
512:     } else if (rblks == 2) {

514:       /* one merge in right branch, right branch done */
515:       /* (note that nothing needs to be done if RBLKS.EQ.1 !) */

517:       ++merged;
518:       startp += lsum;
519:       iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
520:       iwork[icut + (nblks - 1) - merged-1] = start + lblks;
521:       iwork[isize + (nblks - 1) - merged-1] = size - lsum;
522:       iwork[ilsum + (nblks - 1) - merged-1] = ksizes[start + lblks-1];
523:     }
524:     if (rblks <= 0) SETERRQ1(PETSC_COMM_SELF,1,"dibtdc: ERROR in preprocessing - rblks = %d",rblks);

526:     /* continue on the left side */

528:     size = lsum;
529:     blks = lblks;
530:     goto L200;
531:   }

533:   /*  SIZE = IWORK( ISIZE+NBLKS-2 ) */
534:   /*  MAT1 = IWORK( ILSUM+NBLKS-2 ) */

536:   /* Note: after the dimensions SIZE and MAT1 of the last merging */
537:   /* operation have been determined, an upper bound for the workspace */
538:   /* requirements which is independent of how much deflation occurs in */
539:   /* the last merging operation could be determined as follows */
540:   /* (based on (3.15) and (3.19) from UT-CS-00-447): */

542:   /*  IF( MAT1.LE.N/2 ) THEN */
543:   /*     WSPREQ = 3*N + 3/2*( SIZE-MAT1 )**2 + N*N/2 + MAT1*MAT1 */
544:   /*  ELSE */
545:   /*     WSPREQ = 3*N + 3/2*MAT1*MAT1 + N*N/2 + ( SIZE-MAT1 )**2 */
546:   /*  END IF */

548:   /*  IF( LWORK-WSPREQ.LT.0 )THEN */
549:   /*          not enough work space provided */
550:   /*     INFO = -200 - ( WSPREQ-LWORK ) */
551:   /*     RETURN */
552:   /*  END IF */
553:   /*  However, this is not really useful, since the actual check whether */
554:   /*  enough workspace is provided happens in DMERG2.f ! */

556: /* ************************************************************************* */

558:   /* ...Solve subproblems................................... */

560:   /* Divide the matrix into NBLKS submatrices using rank-r */
561:   /* modifications (cuts) and solve for their eigenvalues and */
562:   /* eigenvectors. Initialize index array to sort eigenvalues. */

564:   /* first block: ...................................... */

566:   /*    correction for block 1: D1 - V1 \Sigma1 V1^T */

568:   ksk = ksizes[0];
569:   rp1 = rank[0];

571:   /* initialize the proper part of Z with the diagonal block D1 */
572:   /* (the correction will be made in Z and then the call of DSYEVD will */
573:   /*  overwrite it with the eigenvectors) */

575:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L", &ksk, &ksk, d, &l1d, z, &ldz));

577:   /* copy D1 into WORK (in order to be able to restore it afterwards) */

579:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L", &ksk, &ksk, d, &l1d, work, &ksk));

581:   /* copy V1 into the first RANK(1) columns of D1 and then */
582:   /* multiply with \Sigma1 */

584:   for (i = 0; i < rank[0]; ++i) {
585:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&ksk, &e[(rp1 + i+1)*l1e], &one, &d[i*l1d], &one));
586:     PetscStackCallBLAS("BLASscal",BLASscal_(&ksk, &e[i + rp1*l1e], &d[i*l1d], &one));
587:   }

589:   /* multiply the first RANK( 1 ) columns of D1 with V1^T and */
590:   /* subtract the result from the proper part of Z (previously */
591:   /* initialized with D1) */

593:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N", "T", &ksk, &ksk, rank, &dmone,
594:           d, &l1d, &e[(rank[0]+1)*l1e], &l1e, &done, z, &ldz));

596:   /* restore the original D1 from WORK */

598:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L", &ksk, &ksk, work, &ksk, d, &l1d));

600:   /* eigenanalysis of block 1 (using DSYEVD) */

602:   PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V", "L", &ksk, z, &ldz, ev, work, &lwork, info));
603:   if (*info) SETERRQ1(PETSC_COMM_SELF,1,"dibtdc: Error in DSYEVD for block 1, info = %d",*info);

605:   /* EV( 1: ) contains the eigenvalues in ascending order */
606:   /* (they are returned this way by DSYEVD) */

608:   for (i = 0; i < ksk; ++i) iwork[i] = i+1;

610:     /* intermediate blocks: .............................. */

612:     np = ksk;

614:     /* remaining number of blocks */

616:     if (nblks > 2) {
617:       for (k = 1; k < nblks-1; ++k) {

619:       /* correction for block K: */
620:       /* Dk - U(k-1) \Sigma(k-1) U(k-1)^T - Vk \Sigmak Vk^T */

622:       ksk = ksizes[k];
623:       rp1 = rank[k];

625:       /* initialize the proper part of Z with the diagonal block Dk */
626:       /* (the correction will be made in Z and then the call of DSYEVD will */
627:       /*  overwrite it with the eigenvectors) */

629:       PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L", &ksk, &ksk, &d[k*l1d*l2d], &l1d, &z[np+np*ldz], &ldz));

631:       /* copy Dk into WORK (in order to be able to restore it afterwards) */

633:       PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L", &ksk, &ksk, &d[k*l1d*l2d], &l1d, work, &ksk));

635:       /* copy U(K-1) into the first RANK(K-1) columns of Dk and then */
636:       /* multiply with \Sigma(K-1) */

638:       for (i = 0; i < rank[k-1]; ++i) {
639:         PetscStackCallBLAS("BLAScopy",BLAScopy_(&ksk, &e[(i+(k-1)*l2e)*l1e], &one, &d[(i+k*l2d)*l1d], &one));
640:         PetscStackCallBLAS("BLASscal",BLASscal_(&ksk, &e[i+(rank[k-1]+(k-1)*l2e)*l1e], &d[(i+k*l2d)*l1d], &one));
641:       }

643:       /* multiply the first RANK(K-1) columns of Dk with U(k-1)^T and */
644:       /* subtract the result from the proper part of Z (previously */
645:       /* initialized with Dk) */

647:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N", "T", &ksk, &ksk, &rank[k-1],
648:                     &dmone, &d[k*l1d*l2d],
649:                     &l1d, &e[(k-1)*l1e*l2e], &l1e, &done, &z[np+np*ldz], &ldz));

651:       /* copy Vk into the first RANK(K) columns of Dk and then */
652:       /* multiply with \Sigmak */

654:       for (i = 0; i < rank[k]; ++i) {
655:         PetscStackCallBLAS("BLAScopy",BLAScopy_(&ksk, &e[(rp1+i+1 + k*l2e)*l1e], &one, &d[(i + k*l2d)*l1d], &one));
656:         PetscStackCallBLAS("BLASscal",BLASscal_(&ksk, &e[i + (rp1 + k*l2e)*l1e], &d[(i + k*l2d)*l1d], &one));
657:       }

659:       /* multiply the first RANK(K) columns of Dk with Vk^T and */
660:       /* subtract the result from the proper part of Z (previously */
661:       /* updated with [- U(k-1) \Sigma(k-1) U(k-1)^T] ) */

663:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N", "T", &ksk, &ksk, &rank[k],
664:                     &dmone, &d[k*l1d*l2d], &l1d,
665:                     &e[(rank[k]+1 + k*l2e)*l1e], &l1e, &done, &z[np+np*ldz], &ldz));

667:       /* restore the original Dk from WORK */

669:       PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L", &ksk, &ksk, work, &ksk, &d[k*l1d*l2d], &l1d));

671:       /* eigenanalysis of block K (using dsyevd) */

673:       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V", "L", &ksk, &z[np+np*ldz],
674:                      &ldz, &ev[np], work, &lwork, info));
675:       if (*info) SETERRQ2(PETSC_COMM_SELF,1,"dibtdc: Error in DSYEVD for block %d, info = %d",k,*info);

677:       /* EV( NPP1: ) contains the eigenvalues in ascending order */
678:       /* (they are returned this way by DSYEVD) */

680:       for (i = 0; i < ksk; ++i) iwork[np + i] = i+1;

682:       /* update NP */
683:       np += ksk;
684:     }
685:   }

687:   /* last block: ....................................... */

689:   /*    correction for block NBLKS: */
690:   /*    D(nblks) - U(nblks-1) \Sigma(nblks-1) U(nblks-1)^T */

692:   ksk = ksizes[nblks-1];

694:   /* initialize the proper part of Z with the diagonal block D(nblks) */
695:   /* (the correction will be made in Z and then the call of DSYEVD will */
696:   /* overwrite it with the eigenvectors) */

698:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L", &ksk, &ksk, &d[(nblks-1)*l1d*l2d], &l1d, &z[np+np*ldz], &ldz));

700:   /* copy D(nblks) into WORK (in order to be able to restore it afterwards) */

702:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L", &ksk, &ksk, &d[(nblks-1)*l1d*l2d], &l1d, work, &ksk));

704:   /* copy U(nblks-1) into the first RANK(nblks-1) columns of D(nblks) and then */
705:   /* multiply with \Sigma(nblks-1) */

707:   for (i = 0; i < rank[nblks-2]; ++i) {
708:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&ksk, &e[(i + (nblks-2)*l2e)*l1e],
709:               &one, &d[(i + (nblks-1)*l2d)*l1d], &one));
710:     PetscStackCallBLAS("BLASscal",BLASscal_(&ksk,
711:               &e[i + (rank[nblks-2] + (nblks-2)*l2e)*l1e],
712:               &d[(i + (nblks-1)*l2d)*l1d], &one));
713:   }

715:   /* multiply the first RANK(nblks-1) columns of D(nblks) with U(nblks-1)^T */
716:   /* and subtract the result from the proper part of Z (previously */
717:   /* initialized with D(nblks) ) */

719:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N", "T", &ksk, &ksk, &rank[nblks - 2],
720:           &dmone, &d[(nblks-1)*l1d*l2d], &l1d,
721:           &e[(nblks-2)*l1e*l2e], &l1e, &done, &z[np+np*ldz], &ldz));

723:   /* restore the original D(nblks) from WORK */

725:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("L", &ksk, &ksk, work, &ksk, &d[(nblks-1)*l1d*l2d], &l1d));

727:   /* eigenanalysis of block NBLKS (using dsyevd) */

729:   PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V", "L", &ksk, &z[np+np*ldz], &ldz, &ev[np], work, &lwork, info));
730:   if (*info) SETERRQ2(PETSC_COMM_SELF,1,"dibtdc: Error in DSYEVD for block %d, info = %d",nblks,*info);

732:   /* EV( NPP1: ) contains the eigenvalues in ascending order */
733:   /* (they are returned this way by DSYEVD) */

735:   for (i = 0; i < ksk; ++i) iwork[np + i] = i+1;

737:   /* note that from here on the entire workspace is available again */


740:   /* Perform all the merging operations. */

742:   vstrt = 0;
743:   for (i = 0; i < nblks-1; ++i) {

745:     /* MATSIZ = total size of the current rank RANK modification problem */

747:     matsiz = iwork[isize + i - 1];
748:     np = iwork[istrtp + i - 1];
749:     kbrk = iwork[icut + i - 1];
750:     mat1 = iwork[ilsum + i - 1];
751:     vstrt += np;

753:     for (j = 0; j < rank[kbrk-1]; ++j) {

755:       /* NOTE: The parameter RHO in DMERG2 is modified in DSRTDF */
756:       /*       (multiplied by 2) ! In order not to change the */
757:       /*       singular value stored in E( :, RANK( KBRK )+1, KBRK ), */
758:       /*       we do not pass on this variable as an argument to DMERG2, */
759:       /*       but we assign a separate variable RHO here which is passed */
760:       /*       on to DMERG2. */
761:       /*       Alternative solution in F90: */
762:       /*       pass E( :,RANK( KBRK )+1,KBRK ) to an INTENT( IN ) parameter */
763:       /*       in DMERG2. */

765:       rho = e[j + (rank[kbrk-1] + (kbrk-1)*l2e)*l1e];

767:       /* eigenvectors are accumulated ( JOBZ.EQ.'D' ) */

769:       BDC_dmerg2_(jobz, j+1, matsiz, &ev[np-1], &z[np-1+(np-1)*ldz],
770:                     ldz, &iwork[np-1], &rho, &e[(j + (kbrk-1)*l2e)*l1e],
771:                     ksizes[kbrk], &e[(rank[kbrk-1]+j+1 + (kbrk-1)*l2e)*l1e],
772:                     ksizes[kbrk-1], mat1, work, lwork, &iwork[n], tol, info, 1);
773:                     
774:       if (*info) SETERRQ1(PETSC_COMM_SELF,1,"dibtdc: Error in dmerg2, info = %d",*info);
775:     }

777:     /* at this point all RANK( KBRK ) rank-one modifications corresponding */
778:     /* to the current off-diagonal block are finished. */
779:     /* Move on to the next off-diagonal block. */

781:   }

783:   /* Re-merge the eigenvalues/vectors which were deflated at the final */
784:   /* merging step by sorting all eigenvalues and eigenvectors according */
785:   /* to the permutation stored in IWORK. */

787:   /* copy eigenvalues and eigenvectors in ordered form into WORK */
788:   /* (eigenvalues into WORK( 1:N ), eigenvectors into WORK( N+1:N+1+N^2 ) ) */

790:   for (i = 0; i < n; ++i) {
791:     j = iwork[i];
792:     work[i] = ev[j-1];
793:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&n, &z[(j-1)*ldz], &one, &work[n*(i+1)], &one));
794:   }

796:   /* copy ordered eigenvalues back from WORK( 1:N ) into EV */

798:   PetscStackCallBLAS("BLAScopy",BLAScopy_(&n, work, &one, ev, &one));

800:   /* copy ordered eigenvectors back from WORK( N+1:N+1+N^2 ) into Z */

802:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("A", &n, &n, &work[n], &n, z, &ldz));
803:   return(0);
804: #endif
805: }