Actual source code: hz.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:    HZ iteration for generalized symmetric-indefinite eigenproblem.
 12:    Based on Matlab code from David Watkins.

 14:    References:

 16:        [1] D.S. Watkins, The Matrix Eigenvalue Problem: GR and Krylov Subspace
 17:            Methods, SIAM, 2007.

 19:        [2] M.A. Brebner, J. Grad, "Eigenvalues of Ax = lambda Bx for real
 20:            symmetric matrices A and B computed by reduction to pseudosymmetric
 21:            form and the HR process", Linear Alg. Appl. 43:99-118, 1982.
 22: */

 24: #include <slepc/private/dsimpl.h>
 25: #include <slepcblaslapack.h>

 27: /*
 28:    Sets up a 2-by-2 matrix to eliminate y in the vector [x y]'.
 29:    Transformation is rotator if sygn = 1 and hyperbolic if sygn = -1.
 30: */
 31: static PetscErrorCode UnifiedRotation(PetscReal x,PetscReal y,PetscReal sygn,PetscReal *rot,PetscReal *rcond,PetscBool *swap)
 32: {
 33:   PetscReal nrm,c,s;

 36:   *swap = PETSC_FALSE;
 37:   if (y == 0) {
 38:     rot[0] = 1.0; rot[1] = 0.0; rot[2] = 0.0; rot[3] = 1.0;
 39:     *rcond = 1.0;
 40:   } else {
 41:     nrm = PetscMax(PetscAbs(x),PetscAbs(y));
 42:     c = x/nrm; s = y/nrm;
 43:     if (sygn == 1.0) {  /* set up a rotator */
 44:       nrm = PetscSqrtReal(c*c+s*s);
 45:       c = c/nrm; s = s/nrm;
 46:       /* rot = [c s; -s c]; */
 47:       rot[0] = c; rot[1] = -s; rot[2] = s; rot[3] = c;
 48:       *rcond = 1.0;
 49:     } else if (sygn == -1) {  /* set up a hyperbolic transformation */
 50:       nrm = c*c-s*s;
 51:       if (nrm > 0) nrm = PetscSqrtReal(nrm);
 52:       else if (nrm < 0) {
 53:         nrm = PetscSqrtReal(-nrm);
 54:         *swap = PETSC_TRUE;
 55:       } else SETERRQ(PETSC_COMM_SELF,1,"Breakdown in construction of hyperbolic transformation");
 56:       c = c/nrm; s = s/nrm;
 57:       /* rot = [c -s; -s c]; */
 58:       rot[0] = c; rot[1] = -s; rot[2] = -s; rot[3] = c;
 59:       *rcond = PetscAbs(PetscAbs(s)-PetscAbs(c))/(PetscAbs(s)+PetscAbs(c));
 60:     } else SETERRQ(PETSC_COMM_SELF,1,"Value of sygn sent to transetup must be 1 or -1");
 61:   }
 62:   return(0);
 63: }

 65: static PetscErrorCode HZStep(PetscBLASInt ntop,PetscBLASInt nn,PetscReal tr,PetscReal dt,PetscReal *aa,PetscReal *bb,PetscReal *dd,PetscScalar *uu,PetscInt n,PetscInt ld,PetscBool *flag)
 66: {
 68:   PetscBLASInt   one=1;
 69:   PetscInt       k,jj,ii;
 70:   PetscBLASInt   n_;
 71:   PetscReal      bulge10,bulge20,bulge30,bulge31,bulge41,bulge42;
 72:   PetscReal      sygn,rcond=1.0,worstcond,rot[4],buf[2],t;
 73:   PetscScalar    rtmp;
 74:   PetscBool      swap;

 77:   *flag = PETSC_FALSE;
 78:   worstcond = 1.0;
 79:   PetscBLASIntCast(n,&n_);

 81:   /* Build initial bulge that sets step in motion */
 82:   bulge10 = dd[ntop+1]*(aa[ntop]*(aa[ntop] - dd[ntop]*tr) + dt*dd[ntop]*dd[ntop]) + dd[ntop]*bb[ntop]*bb[ntop];
 83:   bulge20 = bb[ntop]*(dd[ntop+1]*aa[ntop] + dd[ntop]*aa[ntop+1] - dd[ntop]*dd[ntop+1]*tr);
 84:   bulge30 = bb[ntop]*bb[ntop+1]*dd[ntop];
 85:   bulge31 = 0.0;
 86:   bulge41 = 0.0;
 87:   bulge42 = 0.0;

 89:   /* Chase the bulge */
 90:   for (jj=ntop;jj<nn-1;jj++) {

 92:     /* Check for trivial bulge */
 93:     if (jj>ntop && PetscMax(PetscMax(PetscAbs(bulge10),PetscAbs(bulge20)),PetscAbs(bulge30))<PETSC_MACHINE_EPSILON*(PetscAbs(aa[jj]) + PetscAbs(aa[jj+1]))) {
 94:       bb[jj-1] = 0.0;  /* deflate and move on */

 96:     } else { /* carry out the step */

 98:       /* Annihilate tip entry bulge30 */
 99:       if (bulge30 != 0.0) {

101:         /* Make an interchange if necessary to ensure that the
102:            first transformation is othogonal, not hyperbolic.  */
103:         if (dd[jj+1] != dd[jj+2]) { /* make an interchange */
104:           if (dd[jj] != dd[jj+1]) {  /* interchange 1st and 2nd */
105:             buf[0] = bulge20; bulge20 = bulge10; bulge10 = buf[0];
106:             buf[0] = aa[jj]; aa[jj] = aa[jj+1]; aa[jj+1] = buf[0];
107:             buf[0] = bb[jj+1]; bb[jj+1] = bulge31; bulge31 = buf[0];
108:             buf[0] = dd[jj]; dd[jj] = dd[jj+1]; dd[jj+1] = buf[0];
109:             for (k=0;k<n;k++) {
110:               rtmp = uu[k+jj*ld]; uu[k+jj*ld] = uu[k+(jj+1)*ld]; uu[k+(jj+1)*ld] = rtmp;
111:             }
112:           } else {  /* interchange 1st and 3rd */
113:             buf[0] = bulge30; bulge30 = bulge10; bulge10 = buf[0];
114:             buf[0] = aa[jj]; aa[jj] = aa[jj+2]; aa[jj+2] = buf[0];
115:             buf[0] = bb[jj]; bb[jj] = bb[jj+1]; bb[jj+1] = buf[0];
116:             buf[0] = dd[jj]; dd[jj] = dd[jj+2]; dd[jj+2] = buf[0];
117:             if (jj + 2 < nn-1) {
118:               bulge41 = bb[jj+2];
119:               bb[jj+2] = 0;
120:             }
121:             for (k=0;k<n;k++) {
122:               rtmp = uu[k+jj*ld]; uu[k+jj*ld] = uu[k+(jj+2)*ld]; uu[k+(jj+2)*ld] = rtmp;
123:             }
124:           }
125:         }

127:         /* Set up transforming matrix rot. */
128:         UnifiedRotation(bulge20,bulge30,1,rot,&rcond,&swap);

130:         /* Apply transforming matrix rot to T. */
131:         bulge20 = rot[0]*bulge20 + rot[2]*bulge30;
132:         buf[0] = rot[0]*bb[jj] + rot[2]*bulge31;
133:         buf[1] = rot[1]*bb[jj] + rot[3]*bulge31;
134:         bb[jj] = buf[0];
135:         bulge31 = buf[1];
136:         buf[0] = rot[0]*rot[0]*aa[jj+1] + 2.0*rot[0]*rot[2]*bb[jj+1] + rot[2]*rot[2]*aa[jj+2];
137:         buf[1] = rot[1]*rot[1]*aa[jj+1] + 2.0*rot[3]*rot[1]*bb[jj+1] + rot[3]*rot[3]*aa[jj+2];
138:         bb[jj+1] = rot[1]*rot[0]*aa[jj+1] + rot[3]*rot[2]*aa[jj+2] + (rot[3]*rot[0] + rot[1]*rot[2])*bb[jj+1];
139:         aa[jj+1] = buf[0];
140:         aa[jj+2] = buf[1];
141:         if (jj + 2 < nn-1) {
142:           bulge42 = bb[jj+2]*rot[2];
143:           bb[jj+2] = bb[jj+2]*rot[3];
144:         }

146:         /* Accumulate transforming matrix */
147:         PetscStackCallBLAS("BLASrot",BLASMIXEDrot_(&n_,uu+(jj+1)*ld,&one,uu+(jj+2)*ld,&one,&rot[0],&rot[2]));
148:       }

150:       /* Annihilate inner entry bulge20 */
151:       if (bulge20 != 0.0) {

153:         /* Begin by setting up transforming matrix rot */
154:         sygn = dd[jj]*dd[jj+1];
155:         UnifiedRotation(bulge10,bulge20,sygn,rot,&rcond,&swap);
156:         if (rcond<PETSC_MACHINE_EPSILON) {
157:           *flag = PETSC_TRUE;
158:           return(0);
159:         }
160:         if (rcond < worstcond) worstcond = rcond;

162:         /* Apply transforming matrix rot to T */
163:         if (jj > ntop) bb[jj-1] = rot[0]*bulge10 + rot[2]*bulge20;
164:         buf[0] = rot[0]*rot[0]*aa[jj] + 2*rot[0]*rot[2]*bb[jj] + rot[2]*rot[2]*aa[jj+1];
165:         buf[1] = rot[1]*rot[1]*aa[jj] + 2*rot[3]*rot[1]*bb[jj] + rot[3]*rot[3]*aa[jj+1];
166:         bb[jj] = rot[1]*rot[0]*aa[jj] + rot[3]*rot[2]*aa[jj+1] + (rot[3]*rot[0] + rot[1]*rot[2])*bb[jj];
167:         aa[jj] = buf[0];
168:         aa[jj+1] = buf[1];
169:         if (jj + 1 < nn-1) {
170:           /* buf = [ bulge31 bb(jj+1) ] * rot' */
171:           buf[0] = rot[0]*bulge31 + rot[2]*bb[jj+1];
172:           buf[1] = rot[1]*bulge31 + rot[3]*bb[jj+1];
173:           bulge31 = buf[0];
174:           bb[jj+1] = buf[1];
175:         }
176:         if (jj + 2 < nn-1) {
177:           /* buf = [bulge41 bulge42] * rot' */
178:           buf[0] = rot[0]*bulge41 + rot[2]*bulge42;
179:           buf[1] = rot[1]*bulge41 + rot[3]*bulge42;
180:           bulge41 = buf[0];
181:           bulge42 = buf[1];
182:         }

184:         /* Apply transforming matrix rot to D */
185:         if (swap == 1) {
186:           buf[0] = dd[jj]; dd[jj] = dd[jj+1]; dd[jj+1] = buf[0];
187:         }

189:         /* Accumulate transforming matrix, uu(jj:jj+1,:) = rot*uu(jj:jj+1,:) */
190:         if (sygn==1) {
191:           PetscStackCallBLAS("BLASrot",BLASMIXEDrot_(&n_,uu+jj*ld,&one,uu+(jj+1)*ld,&one,&rot[0],&rot[2]));
192:         } else {
193:           if (PetscAbsReal(rot[0])>PetscAbsReal(rot[1])) { /* Type I */
194:             t = rot[1]/rot[0];
195:             for (ii=0;ii<n;ii++) {
196:               uu[jj*ld+ii] = rot[0]*uu[jj*ld+ii] + rot[1]*uu[(jj+1)*ld+ii];
197:               uu[(jj+1)*ld+ii] = t*uu[jj*ld+ii] + uu[(jj+1)*ld+ii]/rot[0];
198:             }
199:           } else { /* Type II */
200:             t = rot[0]/rot[1];
201:             for (ii=0;ii<n;ii++) {
202:               rtmp = uu[jj*ld+ii];
203:               uu[jj*ld+ii] = rot[0]*uu[jj*ld+ii] + rot[1]*uu[(jj+1)*ld+ii];
204:               uu[(jj+1)*ld+ii] = t*uu[jj*ld+ii] + rtmp/rot[1];
205:             }
206:           }
207:         }
208:       }
209:     }

211:     /* Adjust bulge for next step */
212:     bulge10 = bb[jj];
213:     bulge20 = bulge31;
214:     bulge30 = bulge41;
215:     bulge31 = bulge42;
216:     bulge41 = 0.0;
217:     bulge42 = 0.0;
218:   }
219:   return(0);
220: }

222: static PetscErrorCode HZIteration(PetscBLASInt nn,PetscBLASInt cgd,PetscReal *aa,PetscReal *bb,PetscReal *dd,PetscScalar *uu,PetscBLASInt ld)
223: {
225:   PetscBLASInt   j2,one=1;
226:   PetscInt       its,nits,nstop,jj,ntop,nbot,ntry;
227:   PetscReal      htr,det,dis,dif,tn,kt,c,s,tr,dt;
228:   PetscBool      flag=PETSC_FALSE;

231:   its = 0;
232:   nbot = nn-1;
233:   nits = 0;
234:   nstop = 40*(nn - cgd);

236:   while (nbot >= cgd && nits < nstop) {

238:     /* Check for zeros on the subdiagonal */
239:     jj = nbot - 1;
240:     while (jj>=cgd && PetscAbs(bb[jj])>PETSC_MACHINE_EPSILON*(PetscAbs(aa[jj])+PetscAbs(aa[jj+1]))) jj = jj-1;
241:     if (jj>=cgd) bb[jj]=0;
242:     ntop = jj + 1;  /* starting point for step */
243:     if (ntop == nbot) {  /* isolate single eigenvalue */
244:       nbot = ntop - 1;
245:       its = 0;
246:     } else if (ntop+1 == nbot) {  /* isolate pair of eigenvalues */
247:       htr = 0.5*(aa[ntop]*dd[ntop] + aa[nbot]*dd[nbot]);
248:       det = dd[ntop]*dd[nbot]*(aa[ntop]*aa[nbot]-bb[ntop]*bb[ntop]);
249:       dis = htr*htr - det;
250:       if (dis > 0) {  /* distinct real eigenvalues */
251:         if (dd[ntop] == dd[nbot]) {  /* separate the eigenvalues by a Jacobi rotator */
252:           dif = aa[ntop]-aa[nbot];
253:           if (2.0*PetscAbs(bb[ntop])<=dif) {
254:             tn = 2*bb[ntop]/dif;
255:             tn = tn/(1.0 + PetscSqrtScalar(1.0+tn*tn));
256:           } else {
257:             kt = dif/(2.0*bb[ntop]);
258:             tn = PetscSign(kt)/(PetscAbs(kt)+PetscSqrtScalar(1.0+kt*kt));
259:           }
260:           c = 1.0/PetscSqrtScalar(1.0 + tn*tn);
261:           s = c*tn;
262:           aa[ntop] = aa[ntop] + tn*bb[ntop];
263:           aa[nbot] = aa[nbot] - tn*bb[ntop];
264:           bb[ntop] = 0;
265:           j2 = nn-cgd;
266:           PetscStackCallBLAS("BLASrot",BLASMIXEDrot_(&j2,uu+ntop*ld+cgd,&one,uu+nbot*ld+cgd,&one,&c,&s));
267:         }
268:       }
269:       nbot = ntop - 1;
270:     } else {  /* Do an HZ iteration */
271:       its = its + 1;
272:       nits = nits + 1;
273:       tr = aa[nbot-1]*dd[nbot-1] + aa[nbot]*dd[nbot];
274:       dt = dd[nbot-1]*dd[nbot]*(aa[nbot-1]*aa[nbot]-bb[nbot-1]*bb[nbot-1]);
275:       for (ntry=1;ntry<=6;ntry++) {
276:         HZStep(ntop,nbot+1,tr,dt,aa,bb,dd,uu,nn,ld,&flag);
277:         if (!flag) break;
278:         else if (ntry == 6) SETERRQ(PETSC_COMM_SELF,1,"Unable to complete hz step on six tries");
279:         else {
280:           tr = 0.9*tr; dt = 0.81*dt;
281:         }
282:       }
283:     }
284:   }
285:   return(0);
286: }

288: PetscErrorCode DSSolve_GHIEP_HZ(DS ds,PetscScalar *wr,PetscScalar *wi)
289: {
291:   PetscInt       off;
292:   PetscBLASInt   n1,ld;
293:   PetscScalar    *A,*B,*Q;
294:   PetscReal      *d,*e,*s;
295: #if defined(PETSC_USE_COMPLEX)
296:   PetscInt       i;
297: #endif

300: #if !defined(PETSC_USE_COMPLEX)
302: #endif
303:   PetscBLASIntCast(ds->ld,&ld);
304:   n1  = ds->n - ds->l;
305:   off = ds->l + ds->l*ld;
306:   A   = ds->mat[DS_MAT_A];
307:   B   = ds->mat[DS_MAT_B];
308:   Q   = ds->mat[DS_MAT_Q];
309:   d   = ds->rmat[DS_MAT_T];
310:   e   = ds->rmat[DS_MAT_T] + ld;
311:   s   = ds->rmat[DS_MAT_D];
312:   /* Quick return */
313:   if (n1 == 1) {
314:     *(Q+off) = 1;
315:     if (ds->compact) {
316:       wr[ds->l] = d[ds->l]/s[ds->l]; wi[ds->l] = 0.0;
317:     } else {
318:       d[ds->l] = PetscRealPart(A[off]); s[ds->l] = PetscRealPart(B[off]);
319:       wr[ds->l] = d[ds->l]/s[ds->l]; wi[ds->l] = 0.0;
320:     }
321:     return(0);
322:   }
323:   /* Reduce to pseudotriadiagonal form */
324:   DSIntermediate_GHIEP(ds);
325:   HZIteration(ds->n,ds->l,d,e,s,Q,ld);
326:   if (!ds->compact) {
327:     DSSwitchFormat_GHIEP(ds,PETSC_FALSE);
328:   }
329:   /* Undo from diagonal the blocks whith real eigenvalues*/
330:   DSGHIEPRealBlocks(ds);

332:   /* Recover eigenvalues from diagonal */
333:   DSGHIEPComplexEigs(ds,0,ds->n,wr,wi);
334: #if defined(PETSC_USE_COMPLEX)
335:   if (wi) {
336:     for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
337:   }
338: #endif
339:   ds->t = ds->n;
340:   return(0);
341: }