Actual source code: dibtdc.c

slepc-3.12.0 2019-09-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) *info = -1;
112:   else if (n < 3) *info = -2;
113:   else if (blkct < 3) *info = -3;
114:   if (*info == 0) {
115:     ksum = 0;
116:     kchk = 0;
117:     for (i = 0; i < blkct; ++i) {
118:       ksk = bsizes[i];
119:       ksum += ksk;
120:       if (ksk < 1) kchk = 1;
121:     }
122:     if (ksum != n || kchk == 1) *info = -4;
123:   }
124:   if (*info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong argument %d in CUTLR",-(*info));

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

128:   minrnk = n;
129:   for (i = 0; i < blkct-1; ++i) {
130:     if (ranks[i] < minrnk) minrnk = ranks[i];
131:   }

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

135:   nhalf = n / 2;
136:   tmpsum = 0;
137:   mindev = n;
138:   for (i = 0; i < blkct; ++i) {
139:     tmpsum += bsizes[i];
140:     if (ranks[i] == minrnk) {

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

144:       deviat = tmpsum - nhalf;
145:       if (deviat<0) deviat = -deviat;

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

149:       if (deviat < mindev) {
150:         mindev = deviat;
151:         *cut = start + i;
152:         *lblks = i + 1;
153:         *lsum = tmpsum;
154:       }
155:     }
156:   }

158:   if (*cut < start || *cut >= start + blkct - 1) *info = 6;
159:   else if (*lsum < 1 || *lsum >= n) *info = 7;
160:   else if (*lblks < 1 || *lblks >= blkct) *info = 8;
161:   return(0);
162: }

164: PetscErrorCode BDC_dibtdc_(const char *jobz,PetscBLASInt n,PetscBLASInt nblks,
165:         PetscBLASInt *ksizes,PetscReal *d,PetscBLASInt l1d,PetscBLASInt l2d,
166:         PetscReal *e,PetscBLASInt *rank,PetscBLASInt l1e,PetscBLASInt l2e,
167:         PetscReal tol,PetscReal *ev,PetscReal *z,PetscBLASInt ldz,PetscReal *work,
168:         PetscBLASInt lwork,PetscBLASInt *iwork,PetscBLASInt liwork,
169:         PetscBLASInt *info,PetscBLASInt jobz_len)
170: {
171: /*  -- Routine written in LAPACK Version 3.0 style -- */
172: /* *************************************************** */
173: /*     Written by */
174: /*     Michael Moldaschl and Wilfried Gansterer */
175: /*     University of Vienna */
176: /*     last modification: March 16, 2014 */

178: /*     Small adaptations of original code written by */
179: /*     Wilfried Gansterer and Bob Ward, */
180: /*     Department of Computer Science, University of Tennessee */
181: /*     see https://doi.org/10.1137/S1064827501399432 */
182: /* *************************************************** */

184: /*  Purpose */
185: /*  ======= */

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

191: /*  Arguments */
192: /*  ========= */

194: /*  JOBZ    (input) CHARACTER*1 */
195: /*          = 'N':  Compute eigenvalues only (not implemented); */
196: /*          = 'D':  Compute eigenvalues and eigenvectors. */
197: /*                  Eigenvectors are accumulated in the */
198: /*                  divide-and-conquer process. */

200: /*  N      (input) INTEGER */
201: /*         The dimension of the symmetric irreducible block tridiagonal */
202: /*         matrix.  N >= 2. */

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

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

212: /*  D      (input) DOUBLE PRECISION array, dimension (L1D,L2D,NBLKS) */
213: /*         The lower triangular elements of the symmetric diagonal */
214: /*         blocks of the block tridiagonal matrix.  Elements of the top */
215: /*         left diagonal block, which is of dimension KSIZES(1), are */
216: /*         contained in D(*,*,1); the elements of the next diagonal */
217: /*         block, which is of dimension KSIZES(2), are contained in */
218: /*         D(*,*,2); etc. */

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

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

228: /*  E      (input) DOUBLE PRECISION array, dimension (L1E,L2E,NBLKS-1) */
229: /*         Contains the elements of the scalars (singular values) and */
230: /*         vectors (singular vectors) defining the rank RANK subdiagonal */
231: /*         blocks of the matrix. */
232: /*         E(1:RANK(K),RANK(K)+1,K) holds the RANK(K) scalars, */
233: /*         E(:,1:RANK(K),K) holds the RANK(K) column vectors, and */
234: /*         E(:,RANK(K)+2:2*RANK(K)+1,K) holds the row vectors for the K-th */
235: /*         subdiagonal block. */

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

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

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

249: /*  TOL    (input) DOUBLE PRECISION, TOL <= 1.0D-1 */
250: /*         User specified deflation tolerance for the routine DMERG2. */
251: /*         If ( 1.0D-1 >= TOL >= 20*EPS ) then TOL is used as */
252: /*         the deflation tolerance in DSRTDF. */
253: /*         If ( TOL < 20*EPS ) then the standard deflation tolerance from */
254: /*         LAPACK is used as the deflation tolerance in DSRTDF. */

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

260: /*  Z      (input/output) DOUBLE PRECISION array, dimension (LDZ, N) */
261: /*         On entry, Z will be the identity matrix. */
262: /*         On exit, Z contains the eigenvectors of the block tridiagonal */
263: /*         matrix. */

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

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

270: /*  LWORK   (input) INTEGER */
271: /*          The dimension of the array WORK. */
272: /*          In order to guarantee correct results in all cases, */
273: /*          LWORK must be at least ( 2*N**2 + 3*N ). In many cases, */
274: /*          less workspace is required. The absolute minimum required is */
275: /*          ( N**2 + 3*N ). */
276: /*          If the workspace provided is not sufficient, the routine will */
277: /*          return a corresponding error code and report how much workspace */
278: /*          was missing (see INFO). */

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

282: /*  LIWORK  (input) INTEGER */
283: /*          The dimension of the array IWORK. */
284: /*          LIWORK must be at least ( 5*N + 3 + 4*NBLKS - 4 ): */
285: /*                 5*KMAX+3 for DSYEVD, 5*N for ????, */
286: /*                 4*NBLKS-4 for the preprocessing (merging order) */
287: /*          Summarizing, the minimum integer workspace needed is */
288: /*          MAX( 5*N, 5*KMAX + 3 ) + 4*NBLKS - 4 */

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

302: /*  Further Details */
303: /*  =============== */

305: /*  Based on code written by */
306: /*     Wilfried Gansterer and Bob Ward, */
307: /*     Department of Computer Science, University of Tennessee */

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

311: /*  ===================================================================== */

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

326:   *info = 0;

328:   if (*(unsigned char *)jobz != 'N' && *(unsigned char *)jobz != 'D') *info = -1;
329:   else if (n < 2) *info = -2;
330:   else if (nblks < 2 || nblks > n) *info = -3;
331:   if (*info == 0) {
332:     ksum = 0;
333:     kmax = 0;
334:     kchk = 0;
335:     for (k = 0; k < nblks; ++k) {
336:       ksk = ksizes[k];
337:       ksum += ksk;
338:       if (ksk > kmax) kmax = ksk;
339:       if (ksk < 1) kchk = 1;
340:     }
341:     lwmin = n*n + n * 3;
342:     liwmin = PetscMax(n * 5,kmax * 5 + 3) + 4*nblks - 4;
343:     if (ksum != n || kchk == 1) *info = -4;
344:     else if (l1d < PetscMax(3,kmax)) *info = -6;
345:     else if (l2d < PetscMax(3,kmax)) *info = -7;
346:     else if (l1e < PetscMax(3,2*kmax + 1)) *info = -10;
347:     else if (l2e < PetscMax(3,2*kmax + 1)) *info = -11;
348:     else if (tol > .1) *info = -12;
349:     else if (ldz < PetscMax(1,n)) *info = -15;
350:     else if (lwork < lwmin) *info = -17;
351:     else if (liwork < liwmin) *info = -19;
352:   }
353:   if (*info == 0) {
354:     for (k = 0; k < nblks-1; ++k) {
355:       if (rank[k] > PetscMin(ksizes[k],ksizes[k+1]) || rank[k] < 1) *info = -9;
356:     }
357:   }

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

361: /* **************************************************************************** */

363:   /* ...Preprocessing..................................................... */
364:   /*    Determine the optimal order for merging the subblocks and how much */
365:   /*    workspace will be needed for the merging (determined by the last */
366:   /*    merge). Cutpoints for the merging operations are determined and stored */
367:   /*    in reverse chronological order (starting with the final merging */
368:   /*    operation). */

370:   /*    integer workspace requirements for the preprocessing: */
371:   /*         4*(NBLKS-1) for merging history */
372:   /*         at most 3*(NBLKS-1) for stack */

374:   start = 1;
375:   size = n;
376:   blks = nblks;
377:   merged = 0;
378:   k = 0;

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

384:   istck1 = 0;
385:   istck2 = istck1 + nblks;
386:   istck3 = istck2 + nblks;

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

392:   istrtp = n * 5 + 4;
393:   icut = istrtp + nblks - 1;
394:   isize = icut + nblks - 1;
395:   ilsum = isize + nblks - 1;

397: L200:

399:   if (nblks >= 3) {

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

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

409:   } else {
410:     cut = 1;
411:     lsum = ksizes[0];
412:     lblks = 1;
413:   }

415:   ++merged;
416:   startp = 0;
417:   for (i = 0; i < start-1; ++i) startp += ksizes[i];
418:   iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
419:   iwork[icut + (nblks - 1) - merged-1] = cut;
420:   iwork[isize + (nblks - 1) - merged-1] = size;
421:   iwork[ilsum + (nblks - 1) - merged-1] = lsum;

423:   if (lblks == 2) {

425:     /* one merge in left branch, left branch done */
426:     ++merged;
427:     iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
428:     iwork[icut + (nblks - 1) - merged-1] = start;
429:     iwork[isize + (nblks - 1) - merged-1] = lsum;
430:     iwork[ilsum + (nblks - 1) - merged-1] = ksizes[start-1];
431:   }

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

435:     /* left branch done, continue on the right side */
436:     start += lblks;
437:     size -= lsum;
438:     blks -= lblks;

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

442:     if (blks == 2) {

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

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

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

457:       if (k >= 1) {

459:         /* something left on the stack */
460:         start = iwork[istck1 + k-1];
461:         size = iwork[istck2 + k-1];
462:         blks = iwork[istck3 + k-1];
463:         --k;
464:         goto L200;
465:       } else {

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

470:         /* exit preprocessing */

472:       }
473:     } else {

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

477:       goto L200;
478:     }
479:   } else {

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

484:     rblks = blks - lblks;
485:     if (rblks >= 3) {
486:       ++k;
487:       iwork[istck1 + k-1] = cut + 1;
488:       iwork[istck2 + k-1] = size - lsum;
489:       iwork[istck3 + k-1] = rblks;
490:     } else if (rblks == 2) {

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

495:       ++merged;
496:       startp += lsum;
497:       iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
498:       iwork[icut + (nblks - 1) - merged-1] = start + lblks;
499:       iwork[isize + (nblks - 1) - merged-1] = size - lsum;
500:       iwork[ilsum + (nblks - 1) - merged-1] = ksizes[start + lblks-1];
501:     }
502:     if (rblks <= 0) SETERRQ1(PETSC_COMM_SELF,1,"dibtdc: ERROR in preprocessing - rblks = %d",rblks);

504:     /* continue on the left side */

506:     size = lsum;
507:     blks = lblks;
508:     goto L200;
509:   }

511:   /*  SIZE = IWORK( ISIZE+NBLKS-2 ) */
512:   /*  MAT1 = IWORK( ILSUM+NBLKS-2 ) */

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

520:   /*  IF( MAT1.LE.N/2 ) THEN */
521:   /*     WSPREQ = 3*N + 3/2*( SIZE-MAT1 )**2 + N*N/2 + MAT1*MAT1 */
522:   /*  ELSE */
523:   /*     WSPREQ = 3*N + 3/2*MAT1*MAT1 + N*N/2 + ( SIZE-MAT1 )**2 */
524:   /*  END IF */

526:   /*  IF( LWORK-WSPREQ.LT.0 )THEN */
527:   /*          not enough work space provided */
528:   /*     INFO = -200 - ( WSPREQ-LWORK ) */
529:   /*     RETURN */
530:   /*  END IF */
531:   /*  However, this is not really useful, since the actual check whether */
532:   /*  enough workspace is provided happens in DMERG2.f ! */

534: /* ************************************************************************* */

536:   /* ...Solve subproblems................................... */

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

542:   /* first block: ...................................... */

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

546:   ksk = ksizes[0];
547:   rp1 = rank[0];

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

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

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

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

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

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

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

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

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

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

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

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

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

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

588:     /* intermediate blocks: .............................. */

590:     np = ksk;

592:     /* remaining number of blocks */

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

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

600:       ksk = ksizes[k];
601:       rp1 = rank[k];

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

660:       /* update NP */
661:       np += ksk;
662:     }
663:   }

665:   /* last block: ....................................... */

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

670:   ksk = ksizes[nblks-1];

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

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

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

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

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

685:   for (i = 0; i < rank[nblks-2]; ++i) {
686:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&ksk, &e[(i + (nblks-2)*l2e)*l1e],
687:               &one, &d[(i + (nblks-1)*l2d)*l1d], &one));
688:     PetscStackCallBLAS("BLASscal",BLASscal_(&ksk,
689:               &e[i + (rank[nblks-2] + (nblks-2)*l2e)*l1e],
690:               &d[(i + (nblks-1)*l2d)*l1d], &one));
691:   }

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

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

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

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

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

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

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

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

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


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

720:   vstrt = 0;
721:   for (i = 0; i < nblks-1; ++i) {

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

725:     matsiz = iwork[isize + i - 1];
726:     np = iwork[istrtp + i - 1];
727:     kbrk = iwork[icut + i - 1];
728:     mat1 = iwork[ilsum + i - 1];
729:     vstrt += np;

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

733:       /* NOTE: The parameter RHO in DMERG2 is modified in DSRTDF */
734:       /*       (multiplied by 2) ! In order not to change the */
735:       /*       singular value stored in E( :, RANK( KBRK )+1, KBRK ), */
736:       /*       we do not pass on this variable as an argument to DMERG2, */
737:       /*       but we assign a separate variable RHO here which is passed */
738:       /*       on to DMERG2. */
739:       /*       Alternative solution in F90: */
740:       /*       pass E( :,RANK( KBRK )+1,KBRK ) to an INTENT( IN ) parameter */
741:       /*       in DMERG2. */

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

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

747:       BDC_dmerg2_(jobz, j+1, matsiz, &ev[np-1], &z[np-1+(np-1)*ldz],
748:                     ldz, &iwork[np-1], &rho, &e[(j + (kbrk-1)*l2e)*l1e],
749:                     ksizes[kbrk], &e[(rank[kbrk-1]+j+1 + (kbrk-1)*l2e)*l1e],
750:                     ksizes[kbrk-1], mat1, work, lwork, &iwork[n], tol, info, 1);
751:                     
752:       if (*info) SETERRQ1(PETSC_COMM_SELF,1,"dibtdc: Error in dmerg2, info = %d",*info);
753:     }

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

759:   }

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

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

768:   for (i = 0; i < n; ++i) {
769:     j = iwork[i];
770:     work[i] = ev[j-1];
771:     PetscStackCallBLAS("BLAScopy",BLAScopy_(&n, &z[(j-1)*ldz], &one, &work[n*(i+1)], &one));
772:   }

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

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

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

780:   PetscStackCallBLAS("LAPACKlacpy",LAPACKlacpy_("A", &n, &n, &work[n], &n, z, &ldz));
781:   return(0);
782: #endif
783: }