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: PetscErrorCode BDC_dmerg2_(const char *jobz,PetscBLASInt rkct,PetscBLASInt n, 18: PetscReal *ev,PetscReal *q,PetscBLASInt ldq,PetscBLASInt *indxq, 19: PetscReal *rho,PetscReal *u,PetscBLASInt sbrkp1,PetscReal *v, 20: PetscBLASInt sbrk,PetscBLASInt cutpnt,PetscReal *work,PetscBLASInt lwork, 21: PetscBLASInt *iwork,PetscReal tol,PetscBLASInt *info,PetscBLASInt jobz_len) 22: {
23: /* -- Routine written in LAPACK Version 3.0 style -- */
24: /* *************************************************** */
25: /* Written by */
26: /* Michael Moldaschl and Wilfried Gansterer */
27: /* University of Vienna */
28: /* last modification: March 16, 2014 */
30: /* Small adaptations of original code written by */
31: /* Wilfried Gansterer and Bob Ward, */
32: /* Department of Computer Science, University of Tennessee */
33: /* see https://doi.org/10.1137/S1064827501399432 */
34: /* *************************************************** */
36: /* Purpose */
37: /* ======= */
39: /* DMERG2 computes the updated eigensystem of a diagonal matrix after */
40: /* modification by a rank-one symmetric matrix. The diagonal matrix */
41: /* consists of two diagonal submatrices, and the vectors defining the */
42: /* rank-1 matrix similarly have two underlying subvectors each. */
43: /* The dimension of the first subproblem is CUTPNT, the dimension of */
44: /* the second subproblem is N-CUTPNT. */
46: /* T = Q(in) (EV(in) + RHO * Z*Z') Q'(in) = Q(out) * EV(out) * Q'(out) */
48: /* where Z = Q'[V U']', where V is a row vector and U is a column */
49: /* vector with dimensions corresponding to the two underlying */
50: /* subproblems. */
52: /* The eigenvectors of the original matrix are stored in Q, and the */
53: /* eigenvalues in EV. The algorithm consists of three stages: */
55: /* The first stage consists of deflating the size of the problem */
56: /* when there are multiple eigenvalues or if there is a zero in */
57: /* the Z vector. For each such occurrence the dimension of the */
58: /* secular equation problem is reduced by one. This stage is */
59: /* performed by the routine DSRTDF. */
61: /* The second stage consists of calculating the updated */
62: /* eigenvalues. This is done by finding the roots of the secular */
63: /* equation via the routine DLAED4 (as called by DLAED3M). */
64: /* This routine also calculates the eigenvectors of the current */
65: /* problem. */
67: /* If(JOBZ.EQ.'D') then the final stage consists */
68: /* of computing the updated eigenvectors directly using the updated */
69: /* eigenvalues. The eigenvectors for the current problem are multiplied */
70: /* with the eigenvectors from the overall problem. */
72: /* Arguments */
73: /* ========= */
75: /* JOBZ (input) CHARACTER*1 */
76: /* = 'N': Compute eigenvalues only (not implemented); */
77: /* = 'D': Compute eigenvalues and eigenvectors. */
78: /* Eigenvectors are accumulated in the divide-and-conquer */
79: /* process. */
81: /* RKCT (input) INTEGER */
82: /* The number of the rank modification which is accounted for */
83: /* (RKCT >= 1). Required parameter, because the update operation of the */
84: /* modification vector can be performed much more efficiently */
85: /* if RKCT.EQ.1. In that case, the eigenvector matrix is still */
86: /* block-diagonal. For RKCT.GE.2 the eigenvector matrix for the update */
87: /* operation has filled up and is a full matrix. */
89: /* N (input) INTEGER */
90: /* The dimension of the symmetric block tridiagonal matrix. */
91: /* N >= 0. */
93: /* EV (input/output) DOUBLE PRECISION array, dimension (N) */
94: /* On entry, the eigenvalues (=diagonal values) of the */
95: /* rank-1-perturbed matrix. */
96: /* On exit, the eigenvalues of the repaired matrix. */
98: /* Q (input/output) DOUBLE PRECISION array, dimension (LDQ,N) */
99: /* On entry, the eigenvectors of the rank-1-perturbed matrix. */
100: /* On exit, the eigenvectors of the repaired tridiagonal matrix. */
102: /* LDQ (input) INTEGER */
103: /* The leading dimension of the array Q. LDQ >= max(1,N). */
105: /* INDXQ (input/output) INTEGER array, dimension (N) */
106: /* On entry, the permutation which separately sorts the two */
107: /* subproblems in EV into ascending order. */
108: /* On exit, the permutation which will reintegrate the */
109: /* subproblems back into sorted order, */
110: /* i.e. EV(INDXQ(I = 1, N)) will be in ascending order. */
112: /* RHO (input/output) DOUBLE PRECISION */
113: /* The scalar in the rank-1 perturbation. Modified (multiplied */
114: /* by 2) in DSRTDF. */
116: /* U (input) DOUBLE PRECISION array; dimension (SBRKP1), where SBRKP1 */
117: /* is the size of the first (original) block after CUTPNT. */
118: /* The column vector of the rank-1 subdiagonal connecting the */
119: /* two diagonal subproblems. */
120: /* Theoretically, zero entries might have to be appended after U */
121: /* in order to make it have dimension (N-CUTPNT). However, this */
122: /* is not required because it can be accounted for in the */
123: /* matrix-vector product using the argument SBRKP1. */
125: /* SBRKP1 (input) INTEGER */
126: /* Dimension of the relevant (non-zero) part of the vector U. */
127: /* Equal to the size of the first original block after the */
128: /* breakpoint. */
130: /* V (input) DOUBLE PRECISION array; dimension (SBRK), where SBRK */
131: /* is the size of the last (original) block before CUTPNT. */
132: /* The row vector of the rank-1 subdiagonal connecting the two */
133: /* diagonal subproblems. */
134: /* Theoretically, zero entries might have to be inserted in front */
135: /* of V in order to make it have dimension (CUTPNT). However, this */
136: /* is not required because it can be accounted for in the */
137: /* matrix-vector product using the argument SBRK. */
139: /* SBRK (input) INTEGER */
140: /* Dimension of the relevant (non-zero) part of the vector V. */
141: /* Equal to the size of the last original block before the */
142: /* breakpoint. */
144: /* CUTPNT (input) INTEGER */
145: /* The location of the last eigenvalue of the leading diagonal */
146: /* block. min(1,N) <= CUTPNT <= max(1,N). */
148: /* WORK (workspace) DOUBLE PRECISION array, dimension (LWORK) */
150: /* LWORK (input) INTEGER */
151: /* The dimension of the array WORK. */
152: /* In order to guarantee correct results in all cases, */
153: /* LWORK must be at least (2*N**2 + 3*N). In many cases, */
154: /* less workspace is required. The absolute minimum required is */
155: /* (N**2 + 3*N). */
156: /* If the workspace provided is not sufficient, the routine will */
157: /* return a corresponding error code and report how much workspace */
158: /* was missing (see INFO). */
159: /* NOTE: This parameter is needed for determining whether enough */
160: /* workspace is provided, and, if not, for computing how */
161: /* much workspace is needed. */
163: /* IWORK (workspace) INTEGER array, dimension (4*N) */
165: /* TOL (input) DOUBLE PRECISION */
166: /* User specified deflation tolerance for the routine DSRTDF. */
168: /* INFO (output) INTEGER */
169: /* = 0: successful exit. */
170: /* < -200: not enough workspace */
171: /* ABS(INFO + 200) numbers have to be stored in addition */
172: /* to the workspace provided, otherwise some eigenvectors */
173: /* will be incorrect. */
174: /* < 0, > -99: if INFO.EQ.-i, the i-th argument had an */
175: /* illegal value. */
176: /* > 0: if INFO.EQ.1, an eigenvalue did not converge */
177: /* if INFO.EQ.2, the deflation counters DZ and DE do not sum */
178: /* up to the total number N-K of components */
179: /* deflated */
181: /* Further Details */
182: /* =============== */
184: /* Based on code written by */
185: /* Wilfried Gansterer and Bob Ward, */
186: /* Department of Computer Science, University of Tennessee */
188: /* Based on the design of the LAPACK code Dlaed1.f written by Jeff */
189: /* Rutter, Computer Science Division, University of California at */
190: /* Berkeley, and modified by Francoise Tisseur, University of Tennessee. */
192: /* ===================================================================== */
194: PetscBLASInt i, k, n1, n2, de, is, dz, iw, iz, iq2, nmc, cpp1;
195: PetscBLASInt indx, indxc, indxp, lwmin, idlmda;
196: PetscBLASInt spneed, coltyp, tmpcut, i__1, i__2, one=1, mone=-1;
197: char defl[1];
198: PetscReal done = 1.0, dzero = 0.0;
200: *info = 0;
201: lwmin = n*n + n * 3;
203: if (n < 0) *info = -3;
204: else if (ldq < PetscMax(1,n)) *info = -6;
205: else if (cutpnt < PetscMin(1,n) || cutpnt > PetscMax(1,n)) *info = -13;
206: else if (lwork < lwmin) *info = -15;
209: /* **************************************************************************** */
211: /* Quick return if possible */
213: if (n == 0) PetscFunctionReturn(0);
215: /* **************************************************************************** */
217: /* The following values are integer pointers which indicate */
218: /* the portion of the workspace used by a particular array in DSRTDF */
219: /* and DLAED3M. */
221: iz = 0;
222: idlmda = iz + n;
223: iw = idlmda + n;
224: iq2 = iw + n;
225: is = iq2 + n * n;
227: /* After the pointer IS the matrix S is stored and read in WORK */
228: /* in the routine DLAED3M. */
230: indx = 0;
231: indxc = indx + n;
232: coltyp = indxc + n;
233: indxp = coltyp + n;
235: /* If eigenvectors are to be accumulated in the divide-and-conquer */
236: /* process (JOBZ.EQ.'D') form the z-vector which consists of */
237: /* Q_1^T * V and Q_2^T * U. */
239: cpp1 = cutpnt + 1;
240: if (rkct == 1) {
242: /* for the first rank modification the eigenvector matrix has */
243: /* special block-diagonal structure and therefore Q_1^T * V and */
244: /* Q_2^T * U can be formed separately. */
246: PetscStackCallBLAS("BLASgemv",BLASgemv_("T", &sbrk, &cutpnt, &done,
247: &q[cutpnt - sbrk], &ldq, v, &one, &dzero, &work[iz], &one));
248: nmc = n - cutpnt;
249: PetscStackCallBLAS("BLASgemv",BLASgemv_("T", &sbrkp1, &nmc, &done,
250: &q[cpp1-1 + (cpp1-1)*ldq], &ldq, u,
251: &one, &dzero, &work[iz + cutpnt], &one));
253: } else {
255: /* for the higher rank modifications, the vectors V and U */
256: /* have to be multiplied with the full eigenvector matrix */
258: PetscStackCallBLAS("BLASgemv",BLASgemv_("T", &sbrk, &n, &done,
259: &q[cutpnt - sbrk], &ldq, v, &one, &dzero, &work[iz], &one));
260: PetscStackCallBLAS("BLASgemv",BLASgemv_("T", &sbrkp1, &n, &done, &q[cpp1-1],
261: &ldq, u, &one, &done, &work[iz], &one));
263: }
265: /* **************************************************************************** */
267: /* Deflate eigenvalues. */
269: if (rkct == 1) {
271: /* for the first rank modification we need the actual cutpoint */
273: n1 = cutpnt;
274: tmpcut = cutpnt;
275: } else {
277: /* for the later rank modifications there is no actual cutpoint any more */
279: n1 = n;
281: /* The original value of CUTPNT has to be preserved for the next time */
282: /* this subroutine is called (therefore, CUTPNT is an INPUT parameter */
283: /* and not to be changed). Thus, assign N to TMPCUT and use the local */
284: /* variable TMPCUT from now on for the cut point. */
286: tmpcut = n;
287: }
289: /* initialize the flag DEFL (indicates whether deflation occurred - */
290: /* this information is needed later in DLAED3M) */
292: *(unsigned char *)defl = '0';
294: /* call DSRTDF for deflation */
296: PetscCall(BDC_dsrtdf_(&k, n, n1, ev, q, ldq, indxq, rho, &work[iz],
297: &work[idlmda], &work[iw], &work[iq2], &iwork[indx],
298: &iwork[indxc], &iwork[indxp], &iwork[coltyp], tol, &dz, &de, info));
301: if (k < n) {
303: /* ....some deflation occurred in dsrtdf, set the flag DEFL */
304: /* (needed in DLAED3M.f, since Givens rotations need to be */
305: /* applied to the eigenvector matrix only if some deflation */
306: /* happened) */
308: *(unsigned char *)defl = '1';
309: }
311: /* **************************************************************************** */
313: /* Solve the Secular Equation. */
315: if (k != 0) {
317: /* ....not everything was deflated. */
319: /* ....check whether enough workspace is available: */
321: /* Note that the following (upper) bound SPNEED for the workspace */
322: /* requirements should also hold in the extreme case TMPCUT=N, */
323: /* which happens for every rank modification after the first one. */
325: i__1 = (iwork[coltyp] + iwork[coltyp+1]) * k;
326: i__2 = (iwork[coltyp+1] + iwork[coltyp + 2]) * k;
327: spneed = is + PetscMax(i__1,i__2) - 1;
331: /* calling DLAED3M for solving the secular equation. */
333: PetscCall(BDC_dlaed3m_(jobz, defl, k, n, tmpcut, ev, q, ldq,
334: *rho, &work[idlmda], &work[iq2], &iwork[indxc], &iwork[coltyp],
335: &work[iw], &work[is], info, 1, 1));
338: /* Prepare the INDXQ sorting permutation. */
340: n1 = k;
341: n2 = n - k;
342: PetscStackCallBLAS("LAPACKlamrg",LAPACKlamrg_(&n1, &n2, ev, &one, &mone, indxq));
343: if (k == 0) for (i = 0; i < n; ++i) indxq[i] = i+1;
345: } else {
347: /* ....everything was deflated (K.EQ.0) */
349: for (i = 0; i < n; ++i) indxq[i] = i+1;
350: }
351: PetscFunctionReturn(0);
352: }