1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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;
105: *info = 0;
106: *lblks = 1;
107: *lsum = 1;
108: *cut = start;
110: if (start < 1) *info = -1;
111: else if (n < 3) *info = -2;
112: else if (blkct < 3) *info = -3;
113: if (*info == 0) {
114: ksum = 0;
115: kchk = 0;
116: for (i = 0; i < blkct; ++i) {
117: ksk = bsizes[i];
118: ksum += ksk;
119: if (ksk < 1) kchk = 1;
120: }
121: if (ksum != n || kchk == 1) *info = -4;
122: }
125: /* determine smallest rank in the range considered */
127: minrnk = n;
128: for (i = 0; i < blkct-1; ++i) {
129: if (ranks[i] < minrnk) minrnk = ranks[i];
130: }
132: /* determine best cut among those with smallest rank */
134: nhalf = n / 2;
135: tmpsum = 0;
136: mindev = n;
137: for (i = 0; i < blkct; ++i) {
138: tmpsum += bsizes[i];
139: if (ranks[i] == minrnk) {
141: /* determine deviation from "optimal" cut NHALF */
143: deviat = tmpsum - nhalf;
144: if (deviat<0) deviat = -deviat;
146: /* compare to best deviation so far */
148: if (deviat < mindev) {
149: mindev = deviat;
150: *cut = start + i;
151: *lblks = i + 1;
152: *lsum = tmpsum;
153: }
154: }
155: }
157: if (*cut < start || *cut >= start + blkct - 1) *info = 6;
158: else if (*lsum < 1 || *lsum >= n) *info = 7;
159: else if (*lblks < 1 || *lblks >= blkct) *info = 8;
160: PetscFunctionReturn(0);
161: }
163: PetscErrorCode BDC_dibtdc_(const char *jobz,PetscBLASInt n,PetscBLASInt nblks,164: PetscBLASInt *ksizes,PetscReal *d,PetscBLASInt l1d,PetscBLASInt l2d,165: PetscReal *e,PetscBLASInt *rank,PetscBLASInt l1e,PetscBLASInt l2e,166: PetscReal tol,PetscReal *ev,PetscReal *z,PetscBLASInt ldz,PetscReal *work,167: PetscBLASInt lwork,PetscBLASInt *iwork,PetscBLASInt liwork,168: PetscBLASInt *info,PetscBLASInt jobz_len)169: {
170: /* -- Routine written in LAPACK Version 3.0 style -- */
171: /* *************************************************** */
172: /* Written by */
173: /* Michael Moldaschl and Wilfried Gansterer */
174: /* University of Vienna */
175: /* last modification: March 16, 2014 */
177: /* Small adaptations of original code written by */
178: /* Wilfried Gansterer and Bob Ward, */
179: /* Department of Computer Science, University of Tennessee */
180: /* see https://doi.org/10.1137/S1064827501399432 */
181: /* *************************************************** */
183: /* Purpose */
184: /* ======= */
186: /* DIBTDC computes all eigenvalues and corresponding eigenvectors of a */
187: /* symmetric irreducible block tridiagonal matrix with rank RANK matrices */
188: /* as the subdiagonal blocks using a block divide and conquer method. */
190: /* Arguments */
191: /* ========= */
193: /* JOBZ (input) CHARACTER*1 */
194: /* = 'N': Compute eigenvalues only (not implemented); */
195: /* = 'D': Compute eigenvalues and eigenvectors. */
196: /* Eigenvectors are accumulated in the */
197: /* divide-and-conquer process. */
199: /* N (input) INTEGER */
200: /* The dimension of the symmetric irreducible block tridiagonal */
201: /* matrix. N >= 2. */
203: /* NBLKS (input) INTEGER, 2 <= NBLKS <= N */
204: /* The number of diagonal blocks in the matrix. */
206: /* KSIZES (input) INTEGER array, dimension (NBLKS) */
207: /* The dimension of the square diagonal blocks from top left */
208: /* to bottom right. KSIZES(I) >= 1 for all I, and the sum of */
209: /* KSIZES(I) for I = 1 to NBLKS has to be equal to N. */
211: /* D (input) DOUBLE PRECISION array, dimension (L1D,L2D,NBLKS) */
212: /* The lower triangular elements of the symmetric diagonal */
213: /* blocks of the block tridiagonal matrix. Elements of the top */
214: /* left diagonal block, which is of dimension KSIZES(1), are */
215: /* contained in D(*,*,1); the elements of the next diagonal */
216: /* block, which is of dimension KSIZES(2), are contained in */
217: /* D(*,*,2); etc. */
219: /* L1D (input) INTEGER */
220: /* The leading dimension of the array D. L1D >= max(3,KMAX), */
221: /* where KMAX is the dimension of the largest diagonal block. */
223: /* L2D (input) INTEGER */
224: /* The second dimension of the array D. L2D >= max(3,KMAX), */
225: /* where KMAX is as stated in L1D above. */
227: /* E (input) DOUBLE PRECISION array, dimension (L1E,L2E,NBLKS-1) */
228: /* Contains the elements of the scalars (singular values) and */
229: /* vectors (singular vectors) defining the rank RANK subdiagonal */
230: /* blocks of the matrix. */
231: /* E(1:RANK(K),RANK(K)+1,K) holds the RANK(K) scalars, */
232: /* E(:,1:RANK(K),K) holds the RANK(K) column vectors, and */
233: /* E(:,RANK(K)+2:2*RANK(K)+1,K) holds the row vectors for the K-th */
234: /* subdiagonal block. */
236: /* RANK (input) INTEGER array, dimension (NBLKS-1). */
237: /* The ranks of all the subdiagonal blocks contained in the array E. */
238: /* RANK(K) <= MIN(KSIZES(K), KSIZES(K+1)) */
240: /* L1E (input) INTEGER */
241: /* The leading dimension of the array E. L1E >= max(3,2*KMAX+1), */
242: /* where KMAX is as stated in L1D above. */
244: /* L2E (input) INTEGER */
245: /* The second dimension of the array E. L2E >= max(3,2*KMAX+1), */
246: /* where KMAX is as stated in L1D above. */
248: /* TOL (input) DOUBLE PRECISION, TOL <= 1.0D-1 */
249: /* User specified deflation tolerance for the routine DMERG2. */
250: /* If (1.0D-1 >= TOL >= 20*EPS) then TOL is used as */
251: /* the deflation tolerance in DSRTDF. */
252: /* If (TOL < 20*EPS) then the standard deflation tolerance from */
253: /* LAPACK is used as the deflation tolerance in DSRTDF. */
255: /* EV (output) DOUBLE PRECISION array, dimension (N) */
256: /* If INFO = 0, then EV contains the eigenvalues of the */
257: /* symmetric block tridiagonal matrix in ascending order. */
259: /* Z (input/output) DOUBLE PRECISION array, dimension (LDZ, N) */
260: /* On entry, Z will be the identity matrix. */
261: /* On exit, Z contains the eigenvectors of the block tridiagonal */
262: /* matrix. */
264: /* LDZ (input) INTEGER */
265: /* The leading dimension of the array Z. LDZ >= max(1,N). */
267: /* WORK (workspace) DOUBLE PRECISION array, dimension (LWORK) */
269: /* LWORK (input) INTEGER */
270: /* The dimension of the array WORK. */
271: /* In order to guarantee correct results in all cases, */
272: /* LWORK must be at least (2*N**2 + 3*N). In many cases, */
273: /* less workspace is required. The absolute minimum required is */
274: /* (N**2 + 3*N). */
275: /* If the workspace provided is not sufficient, the routine will */
276: /* return a corresponding error code and report how much workspace */
277: /* was missing (see INFO). */
279: /* IWORK (workspace) INTEGER array, dimension (LIWORK) */
281: /* LIWORK (input) INTEGER */
282: /* The dimension of the array IWORK. */
283: /* LIWORK must be at least (5*N + 3 + 4*NBLKS - 4): */
284: /* 5*KMAX+3 for DSYEVD, 5*N for ????, */
285: /* 4*NBLKS-4 for the preprocessing (merging order) */
286: /* Summarizing, the minimum integer workspace needed is */
287: /* MAX(5*N, 5*KMAX + 3) + 4*NBLKS - 4 */
289: /* INFO (output) INTEGER */
290: /* = 0: successful exit. */
291: /* < 0, > -99: illegal arguments. */
292: /* if INFO = -i, the i-th argument had an illegal value. */
293: /* = -99: error in the preprocessing (call of CUTLR). */
294: /* < -200: not enough workspace. Space for ABS(INFO + 200) */
295: /* numbers is required in addition to the workspace provided, */
296: /* otherwise some eigenvectors will be incorrect. */
297: /* > 0: The algorithm failed to compute an eigenvalue while */
298: /* working on the submatrix lying in rows and columns */
299: /* INFO/(N+1) through mod(INFO,N+1). */
301: /* Further Details */
302: /* =============== */
304: /* Based on code written by */
305: /* Wilfried Gansterer and Bob Ward, */
306: /* Department of Computer Science, University of Tennessee */
308: /* This routine is comparable to Dlaed0.f from LAPACK. */
310: /* ===================================================================== */
312: PetscBLASInt i, j, k, np, rp1, ksk, one=1;
313: PetscBLASInt cut, mat1, kchk, kbrk, blks, kmax, icut, size, ksum, lsum;
314: PetscBLASInt lblks, rblks, isize, lwmin, ilsum;
315: PetscBLASInt start, istck1, istck2, istck3, merged;
316: PetscBLASInt liwmin, matsiz, startp, istrtp;
317: PetscReal rho, done=1.0, dmone=-1.0;
319: *info = 0;
321: if (*(unsigned char *)jobz != 'N' && *(unsigned char *)jobz != 'D') *info = -1;
322: else if (n < 2) *info = -2;
323: else if (nblks < 2 || nblks > n) *info = -3;
324: if (*info == 0) {
325: ksum = 0;
326: kmax = 0;
327: kchk = 0;
328: for (k = 0; k < nblks; ++k) {
329: ksk = ksizes[k];
330: ksum += ksk;
331: if (ksk > kmax) kmax = ksk;
332: if (ksk < 1) kchk = 1;
333: }
334: lwmin = n*n + n * 3;
335: liwmin = PetscMax(n * 5,kmax * 5 + 3) + 4*nblks - 4;
336: if (ksum != n || kchk == 1) *info = -4;
337: else if (l1d < PetscMax(3,kmax)) *info = -6;
338: else if (l2d < PetscMax(3,kmax)) *info = -7;
339: else if (l1e < PetscMax(3,2*kmax + 1)) *info = -10;
340: else if (l2e < PetscMax(3,2*kmax + 1)) *info = -11;
341: else if (tol > .1) *info = -12;
342: else if (ldz < PetscMax(1,n)) *info = -15;
343: else if (lwork < lwmin) *info = -17;
344: else if (liwork < liwmin) *info = -19;
345: }
346: if (*info == 0) {
347: for (k = 0; k < nblks-1; ++k) {
348: if (rank[k] > PetscMin(ksizes[k],ksizes[k+1]) || rank[k] < 1) *info = -9;
349: }
350: }
354: /* **************************************************************************** */
356: /* ...Preprocessing..................................................... */
357: /* Determine the optimal order for merging the subblocks and how much */
358: /* workspace will be needed for the merging (determined by the last */
359: /* merge). Cutpoints for the merging operations are determined and stored */
360: /* in reverse chronological order (starting with the final merging */
361: /* operation). */
363: /* integer workspace requirements for the preprocessing: */
364: /* 4*(NBLKS-1) for merging history */
365: /* at most 3*(NBLKS-1) for stack */
367: start = 1;
368: size = n;
369: blks = nblks;
370: merged = 0;
371: k = 0;
373: /* integer workspace used for the stack is not needed any more after the */
374: /* preprocessing and therefore can use part of the 5*N */
375: /* integer workspace needed later on in the code */
377: istck1 = 0;
378: istck2 = istck1 + nblks;
379: istck3 = istck2 + nblks;
381: /* integer workspace used for storing the order of merges starts AFTER */
382: /* the integer workspace 5*N+3 which is needed later on in the code */
383: /* (5*KMAX+3 for DSYEVD, 4*N in DMERG2) */
385: istrtp = n * 5 + 4;
386: icut = istrtp + nblks - 1;
387: isize = icut + nblks - 1;
388: ilsum = isize + nblks - 1;
390: L200:392: if (nblks >= 3) {
394: /* Determine the cut point. Note that in the routine CUTLR it is */
395: /* chosen such that it yields the best balanced merging operation */
396: /* among all the rank modifications with minimum rank. */
398: cutlr_(start, size, blks, &ksizes[start-1], &rank[start-1], &cut, &lsum, &lblks, info);
401: } else {
402: cut = 1;
403: lsum = ksizes[0];
404: lblks = 1;
405: }
407: ++merged;
408: startp = 0;
409: for (i = 0; i < start-1; ++i) startp += ksizes[i];
410: iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
411: iwork[icut + (nblks - 1) - merged-1] = cut;
412: iwork[isize + (nblks - 1) - merged-1] = size;
413: iwork[ilsum + (nblks - 1) - merged-1] = lsum;
415: if (lblks == 2) {
417: /* one merge in left branch, left branch done */
418: ++merged;
419: iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
420: iwork[icut + (nblks - 1) - merged-1] = start;
421: iwork[isize + (nblks - 1) - merged-1] = lsum;
422: iwork[ilsum + (nblks - 1) - merged-1] = ksizes[start-1];
423: }
425: if (lblks == 1 || lblks == 2) {
427: /* left branch done, continue on the right side */
428: start += lblks;
429: size -= lsum;
430: blks -= lblks;
434: if (blks == 2) {
436: /* one merge in right branch, right branch done */
437: ++merged;
438: startp += lsum;
439: iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
440: iwork[icut + (nblks - 1) - merged-1] = start;
441: iwork[isize + (nblks - 1) - merged-1] = size;
442: iwork[ilsum + (nblks - 1) - merged-1] = ksizes[start-1];
443: }
445: if (blks == 1 || blks == 2) {
447: /* get the next subproblem from the stack or finished */
449: if (k >= 1) {
451: /* something left on the stack */
452: start = iwork[istck1 + k-1];
453: size = iwork[istck2 + k-1];
454: blks = iwork[istck3 + k-1];
455: --k;
456: goto L200;
457: } else {
459: /* nothing left on the stack */
462: /* exit preprocessing */
464: }
465: } else {
467: /* BLKS.GE.3, and therefore analyze the right side */
469: goto L200;
470: }
471: } else {
473: /* LBLKS.GE.3, and therefore check the right side and */
474: /* put it on the stack if required */
476: rblks = blks - lblks;
477: if (rblks >= 3) {
478: ++k;
479: iwork[istck1 + k-1] = cut + 1;
480: iwork[istck2 + k-1] = size - lsum;
481: iwork[istck3 + k-1] = rblks;
482: } else if (rblks == 2) {
484: /* one merge in right branch, right branch done */
485: /* (note that nothing needs to be done if RBLKS.EQ.1 !) */
487: ++merged;
488: startp += lsum;
489: iwork[istrtp + (nblks - 1) - merged-1] = startp + 1;
490: iwork[icut + (nblks - 1) - merged-1] = start + lblks;
491: iwork[isize + (nblks - 1) - merged-1] = size - lsum;
492: iwork[ilsum + (nblks - 1) - merged-1] = ksizes[start + lblks-1];
493: }
496: /* continue on the left side */
498: size = lsum;
499: blks = lblks;
500: goto L200;
501: }
503: /* SIZE = IWORK(ISIZE+NBLKS-2) */
504: /* MAT1 = IWORK(ILSUM+NBLKS-2) */
506: /* Note: after the dimensions SIZE and MAT1 of the last merging */
507: /* operation have been determined, an upper bound for the workspace */
508: /* requirements which is independent of how much deflation occurs in */
509: /* the last merging operation could be determined as follows */
510: /* (based on (3.15) and (3.19) from UT-CS-00-447): */
512: /* IF(MAT1.LE.N/2) THEN */
513: /* WSPREQ = 3*N + 3/2*(SIZE-MAT1)**2 + N*N/2 + MAT1*MAT1 */
514: /* ELSE */
515: /* WSPREQ = 3*N + 3/2*MAT1*MAT1 + N*N/2 + (SIZE-MAT1)**2 */
516: /* END IF */
518: /* IF(LWORK-WSPREQ.LT.0)THEN */
519: /* not enough work space provided */
520: /* INFO = -200 - (WSPREQ-LWORK) */
521: /* RETURN */
522: /* END IF */
523: /* However, this is not really useful, since the actual check whether */
524: /* enough workspace is provided happens in DMERG2.f ! */
526: /* ************************************************************************* */
528: /* ...Solve subproblems................................... */
530: /* Divide the matrix into NBLKS submatrices using rank-r */
531: /* modifications (cuts) and solve for their eigenvalues and */
532: /* eigenvectors. Initialize index array to sort eigenvalues. */
534: /* first block: ...................................... */
536: /* correction for block 1: D1 - V1 \Sigma1 V1^T */
538: ksk = ksizes[0];
539: rp1 = rank[0];
541: /* initialize the proper part of Z with the diagonal block D1 */
542: /* (the correction will be made in Z and then the call of DSYEVD will */
543: /* overwrite it with the eigenvectors) */
545: for (j=0;j<ksk;j++) for (i=j;i<ksk;i++) z[i+j*ldz] = d[i+j*l1d];
547: /* copy D1 into WORK (in order to be able to restore it afterwards) */
549: for (j=0;j<ksk;j++) for (i=j;i<ksk;i++) work[i+j*ksk] = d[i+j*l1d];
551: /* copy V1 into the first RANK(1) columns of D1 and then */
552: /* multiply with \Sigma1 */
554: for (i = 0; i < rank[0]; ++i) {
555: PetscStackCallBLAS("BLAScopy",BLAScopy_(&ksk, &e[(rp1 + i+1)*l1e], &one, &d[i*l1d], &one));
556: PetscStackCallBLAS("BLASscal",BLASscal_(&ksk, &e[i + rp1*l1e], &d[i*l1d], &one));
557: }
559: /* multiply the first RANK(1) columns of D1 with V1^T and */
560: /* subtract the result from the proper part of Z (previously */
561: /* initialized with D1) */
563: PetscStackCallBLAS("BLASgemm",BLASgemm_("N", "T", &ksk, &ksk, rank, &dmone,
564: d, &l1d, &e[(rank[0]+1)*l1e], &l1e, &done, z, &ldz));
566: /* restore the original D1 from WORK */
568: for (j=0;j<ksk;j++) for (i=j;i<ksk;i++) d[i+j*l1d] = work[i+j*ksk];
570: /* eigenanalysis of block 1 (using DSYEVD) */
572: PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V", "L", &ksk, z, &ldz, ev, work, &lwork, info));
573: SlepcCheckLapackInfo("syev",*info);
575: /* EV(1:) contains the eigenvalues in ascending order */
576: /* (they are returned this way by DSYEVD) */
578: for (i = 0; i < ksk; ++i) iwork[i] = i+1;
580: /* intermediate blocks: .............................. */
582: np = ksk;
584: /* remaining number of blocks */
586: if (nblks > 2) {
587: for (k = 1; k < nblks-1; ++k) {
589: /* correction for block K: */
590: /* Dk - U(k-1) \Sigma(k-1) U(k-1)^T - Vk \Sigmak Vk^T */
592: ksk = ksizes[k];
593: rp1 = rank[k];
595: /* initialize the proper part of Z with the diagonal block Dk */
596: /* (the correction will be made in Z and then the call of DSYEVD will */
597: /* overwrite it with the eigenvectors) */
599: for (j=0;j<ksk;j++) for (i=j;i<ksk;i++) z[np+i+(np+j)*ldz] = d[i+(j+k*l2d)*l1d];
601: /* copy Dk into WORK (in order to be able to restore it afterwards) */
603: for (j=0;j<ksk;j++) for (i=j;i<ksk;i++) work[i+j*ksk] = d[i+(j+k*l2d)*l1d];
605: /* copy U(K-1) into the first RANK(K-1) columns of Dk and then */
606: /* multiply with \Sigma(K-1) */
608: for (i = 0; i < rank[k-1]; ++i) {
609: PetscStackCallBLAS("BLAScopy",BLAScopy_(&ksk, &e[(i+(k-1)*l2e)*l1e], &one, &d[(i+k*l2d)*l1d], &one));
610: PetscStackCallBLAS("BLASscal",BLASscal_(&ksk, &e[i+(rank[k-1]+(k-1)*l2e)*l1e], &d[(i+k*l2d)*l1d], &one));
611: }
613: /* multiply the first RANK(K-1) columns of Dk with U(k-1)^T and */
614: /* subtract the result from the proper part of Z (previously */
615: /* initialized with Dk) */
617: PetscStackCallBLAS("BLASgemm",BLASgemm_("N", "T", &ksk, &ksk, &rank[k-1],
618: &dmone, &d[k*l1d*l2d],
619: &l1d, &e[(k-1)*l1e*l2e], &l1e, &done, &z[np+np*ldz], &ldz));
621: /* copy Vk into the first RANK(K) columns of Dk and then */
622: /* multiply with \Sigmak */
624: for (i = 0; i < rank[k]; ++i) {
625: PetscStackCallBLAS("BLAScopy",BLAScopy_(&ksk, &e[(rp1+i+1 + k*l2e)*l1e], &one, &d[(i + k*l2d)*l1d], &one));
626: PetscStackCallBLAS("BLASscal",BLASscal_(&ksk, &e[i + (rp1 + k*l2e)*l1e], &d[(i + k*l2d)*l1d], &one));
627: }
629: /* multiply the first RANK(K) columns of Dk with Vk^T and */
630: /* subtract the result from the proper part of Z (previously */
631: /* updated with [- U(k-1) \Sigma(k-1) U(k-1)^T]) */
633: PetscStackCallBLAS("BLASgemm",BLASgemm_("N", "T", &ksk, &ksk, &rank[k],
634: &dmone, &d[k*l1d*l2d], &l1d,
635: &e[(rank[k]+1 + k*l2e)*l1e], &l1e, &done, &z[np+np*ldz], &ldz));
637: /* restore the original Dk from WORK */
639: for (j=0;j<ksk;j++) for (i=j;i<ksk;i++) d[i+(j+k*l2d)*l1d] = work[i+j*ksk];
641: /* eigenanalysis of block K (using dsyevd) */
643: PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V", "L", &ksk, &z[np+np*ldz],
644: &ldz, &ev[np], work, &lwork, info));
645: SlepcCheckLapackInfo("syev",*info);
647: /* EV(NPP1:) contains the eigenvalues in ascending order */
648: /* (they are returned this way by DSYEVD) */
650: for (i = 0; i < ksk; ++i) iwork[np + i] = i+1;
652: /* update NP */
653: np += ksk;
654: }
655: }
657: /* last block: ....................................... */
659: /* correction for block NBLKS: */
660: /* D(nblks) - U(nblks-1) \Sigma(nblks-1) U(nblks-1)^T */
662: ksk = ksizes[nblks-1];
664: /* initialize the proper part of Z with the diagonal block D(nblks) */
665: /* (the correction will be made in Z and then the call of DSYEVD will */
666: /* overwrite it with the eigenvectors) */
668: for (j=0;j<ksk;j++) for (i=j;i<ksk;i++) z[np+i+(np+j)*ldz] = d[i+(j+(nblks-1)*l2d)*l1d];
670: /* copy D(nblks) into WORK (in order to be able to restore it afterwards) */
672: for (j=0;j<ksk;j++) for (i=j;i<ksk;i++) work[i+j*ksk] = d[i+(j+(nblks-1)*l2d)*l1d];
674: /* copy U(nblks-1) into the first RANK(nblks-1) columns of D(nblks) and then */
675: /* multiply with \Sigma(nblks-1) */
677: for (i = 0; i < rank[nblks-2]; ++i) {
678: PetscStackCallBLAS("BLAScopy",BLAScopy_(&ksk, &e[(i + (nblks-2)*l2e)*l1e],
679: &one, &d[(i + (nblks-1)*l2d)*l1d], &one));
680: PetscStackCallBLAS("BLASscal",BLASscal_(&ksk,
681: &e[i + (rank[nblks-2] + (nblks-2)*l2e)*l1e],
682: &d[(i + (nblks-1)*l2d)*l1d], &one));
683: }
685: /* multiply the first RANK(nblks-1) columns of D(nblks) with U(nblks-1)^T */
686: /* and subtract the result from the proper part of Z (previously */
687: /* initialized with D(nblks)) */
689: PetscStackCallBLAS("BLASgemm",BLASgemm_("N", "T", &ksk, &ksk, &rank[nblks - 2],
690: &dmone, &d[(nblks-1)*l1d*l2d], &l1d,
691: &e[(nblks-2)*l1e*l2e], &l1e, &done, &z[np+np*ldz], &ldz));
693: /* restore the original D(nblks) from WORK */
695: for (j=0;j<ksk;j++) for (i=j;i<ksk;i++) d[i+(j+(nblks-1)*l2d)*l1d] = work[i+j*ksk];
697: /* eigenanalysis of block NBLKS (using dsyevd) */
699: PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V", "L", &ksk, &z[np+np*ldz], &ldz, &ev[np], work, &lwork, info));
700: SlepcCheckLapackInfo("syev",*info);
702: /* EV(NPP1:) contains the eigenvalues in ascending order */
703: /* (they are returned this way by DSYEVD) */
705: for (i = 0; i < ksk; ++i) iwork[np + i] = i+1;
707: /* note that from here on the entire workspace is available again */
709: /* Perform all the merging operations. */
711: for (i = 0; i < nblks-1; ++i) {
713: /* MATSIZ = total size of the current rank RANK modification problem */
715: matsiz = iwork[isize + i - 1];
716: np = iwork[istrtp + i - 1];
717: kbrk = iwork[icut + i - 1];
718: mat1 = iwork[ilsum + i - 1];
720: for (j = 0; j < rank[kbrk-1]; ++j) {
722: /* NOTE: The parameter RHO in DMERG2 is modified in DSRTDF */
723: /* (multiplied by 2) ! In order not to change the */
724: /* singular value stored in E(:, RANK(KBRK)+1, KBRK), */
725: /* we do not pass on this variable as an argument to DMERG2, */
726: /* but we assign a separate variable RHO here which is passed */
727: /* on to DMERG2. */
728: /* Alternative solution in F90: */
729: /* pass E(:,RANK(KBRK)+1,KBRK) to an INTENT(IN) parameter */
730: /* in DMERG2. */
732: rho = e[j + (rank[kbrk-1] + (kbrk-1)*l2e)*l1e];
734: /* eigenvectors are accumulated (JOBZ.EQ.'D') */
736: PetscCall(BDC_dmerg2_(jobz, j+1, matsiz, &ev[np-1], &z[np-1+(np-1)*ldz],
737: ldz, &iwork[np-1], &rho, &e[(j + (kbrk-1)*l2e)*l1e],
738: ksizes[kbrk], &e[(rank[kbrk-1]+j+1 + (kbrk-1)*l2e)*l1e],
739: ksizes[kbrk-1], mat1, work, lwork, &iwork[n], tol, info, 1));
741: }
743: /* at this point all RANK(KBRK) rank-one modifications corresponding */
744: /* to the current off-diagonal block are finished. */
745: /* Move on to the next off-diagonal block. */
747: }
749: /* Re-merge the eigenvalues/vectors which were deflated at the final */
750: /* merging step by sorting all eigenvalues and eigenvectors according */
751: /* to the permutation stored in IWORK. */
753: /* copy eigenvalues and eigenvectors in ordered form into WORK */
754: /* (eigenvalues into WORK(1:N), eigenvectors into WORK(N+1:N+1+N^2)) */
756: for (i = 0; i < n; ++i) {
757: j = iwork[i];
758: work[i] = ev[j-1];
759: PetscStackCallBLAS("BLAScopy",BLAScopy_(&n, &z[(j-1)*ldz], &one, &work[n*(i+1)], &one));
760: }
762: /* copy ordered eigenvalues back from WORK(1:N) into EV */
764: PetscStackCallBLAS("BLAScopy",BLAScopy_(&n, work, &one, ev, &one));
766: /* copy ordered eigenvectors back from WORK(N+1:N+1+N^2) into Z */
768: for (j=0;j<n;j++) for (i=0;i<n;i++) z[i+j*ldz] = work[i+(j+1)*n];
769: PetscFunctionReturn(0);
770: }