Actual source code: invit.c

slepc-3.7.4 2017-05-17
Report Typos and Errors
  1: /*

  3:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  4:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  5:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  7:    This file is part of SLEPc.

  9:    SLEPc is free software: you can redistribute it and/or modify it under  the
 10:    terms of version 3 of the GNU Lesser General Public License as published by
 11:    the Free Software Foundation.

 13:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 14:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 15:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 16:    more details.

 18:    You  should have received a copy of the GNU Lesser General  Public  License
 19:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 20:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 21: */
 22: #include <slepc/private/dsimpl.h>
 23: #include <slepcblaslapack.h>

 25: struct HRtr
 26: {
 27:   PetscScalar *data;
 28:   PetscInt    m;
 29:   PetscInt    idx[2];
 30:   PetscInt    n[2];
 31:   PetscScalar tau[2];
 32:   PetscReal   alpha;
 33:   PetscReal   cs;
 34:   PetscReal   sn;
 35:   PetscInt    type;
 36: };

 40: /*
 41:   Generates a hyperbolic rotation
 42:     if x1*x1 - x2*x2 != 0
 43:       r = sqrt(|x1*x1 - x2*x2|)
 44:       c = x1/r  s = x2/r

 46:       | c -s||x1|   |d*r|
 47:       |-s  c||x2| = | 0 |
 48:       where d = 1 for type==1 and -1 for type==2
 49:   Returns the condition number of the reduction
 50: */
 51: static PetscErrorCode HRGen(PetscReal x1,PetscReal x2,PetscInt *type,PetscReal *c,PetscReal *s,PetscReal *r,PetscReal *cond)
 52: {
 53:   PetscReal t,n2,xa,xb;
 54:   PetscInt  type_;

 57:   if (x2==0.0) {
 58:     *r = PetscAbsReal(x1);
 59:     *c = (x1>=0)?1.0:-1.0;
 60:     *s = 0.0;
 61:     if (type) *type = 1;
 62:     return(0);
 63:   }
 64:   if (PetscAbsReal(x1) == PetscAbsReal(x2)) {
 65:     /* hyperbolic rotation doesn't exist */
 66:     *c = 0.0;
 67:     *s = 0.0;
 68:     *r = 0.0;
 69:     if (type) *type = 0;
 70:     *cond = PETSC_MAX_REAL;
 71:     return(0);
 72:   }

 74:   if (PetscAbsReal(x1)>PetscAbsReal(x2)) {
 75:     xa = x1; xb = x2; type_ = 1;
 76:   } else {
 77:     xa = x2; xb = x1; type_ = 2;
 78:   }
 79:   t = xb/xa;
 80:   n2 = PetscAbsReal(1 - t*t);
 81:   *r = PetscSqrtReal(n2)*PetscAbsReal(xa);
 82:   *c = x1/(*r);
 83:   *s = x2/(*r);
 84:   if (type_ == 2) *r *= -1;
 85:   if (type) *type = type_;
 86:   if (cond) *cond = (PetscAbsReal(*c) + PetscAbsReal(*s))/PetscAbsReal(PetscAbsReal(*c) - PetscAbsReal(*s));
 87:   return(0);
 88: }

 92: /*
 93:                                 |c  s|
 94:   Applies an hyperbolic rotator |s  c|
 95:            |c  s|
 96:     [x1 x2]|s  c|
 97: */
 98: static PetscErrorCode HRApply(PetscInt n,PetscScalar *x1,PetscInt inc1,PetscScalar *x2,PetscInt inc2,PetscReal c,PetscReal s)
 99: {
100:   PetscInt    i;
101:   PetscReal   t;
102:   PetscScalar tmp;

105:   if (PetscAbsReal(c)>PetscAbsReal(s)) { /* Type I */
106:     t = s/c;
107:     for (i=0;i<n;i++) {
108:       x1[i*inc1] = c*x1[i*inc1] + s*x2[i*inc2];
109:       x2[i*inc2] = t*x1[i*inc1] + x2[i*inc2]/c;
110:     }
111:   } else { /* Type II */
112:     t = c/s;
113:     for (i=0;i<n;i++) {
114:       tmp = x1[i*inc1];
115:       x1[i*inc1] = c*x1[i*inc1] + s*x2[i*inc2];
116:       x2[i*inc2] = t*x1[i*inc1] + tmp/s;
117:     }
118:   }
119:   return(0);
120: }

124: /*
125:   Reduction to tridiagonal-diagonal form (see F. Tisseur, SIMAX 26(1), 2004).

127:   Input:
128:     A symmetric (only lower triangular part is referred)
129:     s vector +1 and -1 (signature matrix)
130:   Output:
131:     d,e
132:     s
133:     Q s-orthogonal matrix with Q^T*A*Q = T (symmetric tridiagonal matrix)
134: */
135: static PetscErrorCode TridiagDiag_HHR(PetscInt n,PetscScalar *A,PetscInt lda,PetscReal *s,PetscScalar* Q,PetscInt ldq,PetscBool flip,PetscReal *d,PetscReal *e,PetscInt *perm_,PetscScalar *work,PetscReal *rwork,PetscBLASInt *iwork)
136: {
137: #if defined(SLEPC_MISSING_LAPACK_LARFG) || defined(SLEPC_MISSING_LAPACK_LARF)
139:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LARFG/LARF - Lapack routines are unavailable");
140: #else
142:   PetscInt       i,j,k,*ii,*jj,i0=0,ik=0,tmp,type;
143:   PetscInt       nwu=0;
144:   PetscReal      *ss,cond=1.0,cs,sn,r;
145:   PetscScalar    tau,t,*AA;
146:   PetscBLASInt   n0,n1,ni,inc=1,m,n_,lda_,ldq_,*perm;
147:   PetscBool      breakdown = PETSC_TRUE;

150:   if (n<3) {
151:     if (n==1) Q[0]=1;
152:     if (n==2) {
153:       Q[0] = Q[1+ldq] = 1;
154:       Q[1] = Q[ldq] = 0;
155:     }
156:     return(0);
157:   }
158:   PetscBLASIntCast(lda,&lda_);
159:   PetscBLASIntCast(n,&n_);
160:   PetscBLASIntCast(ldq,&ldq_);
161:   ss = rwork;
162:   perm = iwork;
163:   AA = work;
164:   for (i=0;i<n;i++) {
165:     PetscMemcpy(AA+i*n,A+i*lda,n*sizeof(PetscScalar));
166:   }
167:   nwu += n*n;
168:   k=0;
169:   while (breakdown && k<n) {
170:     breakdown = PETSC_FALSE;
171:     /* Classify (and flip) A and s according to sign */
172:     if (flip) {
173:       for (i=0;i<n;i++) {
174:         perm[i] = n-1-perm_[i];
175:         if (perm[i]==0) i0 = i;
176:         if (perm[i]==k) ik = i;
177:       }
178:     } else {
179:       for (i=0;i<n;i++) {
180:         perm[i] = perm_[i];
181:         if (perm[i]==0) i0 = i;
182:         if (perm[i]==k) ik = i;
183:       }
184:     }
185:     perm[ik] = 0;
186:     perm[i0] = k;
187:     i=1;
188:     while (i<n-1 && s[perm[i-1]]==s[perm[0]]) {
189:       if (s[perm[i]]!=s[perm[0]]) {
190:         j=i+1;
191:         while (j<n-1 && s[perm[j]]!=s[perm[0]])j++;
192:         tmp = perm[i]; perm[i] = perm[j]; perm[j] = tmp;
193:       }
194:       i++;
195:     }
196:     for (i=0;i<n;i++) {
197:       ss[i] = s[perm[i]];
198:     }
199:     if (flip) {
200:       ii = &j;
201:       jj = &i;
202:     } else {
203:       ii = &i;
204:       jj = &j;
205:     }
206:     for (i=0;i<n;i++)
207:       for (j=0;j<n;j++)
208:         A[i+j*lda] = AA[perm[*ii]+perm[*jj]*n];
209:     /* Initialize Q */
210:     for (i=0;i<n;i++) {
211:       PetscMemzero(Q+i*ldq,n*sizeof(PetscScalar));
212:       Q[perm[i]+i*ldq] = 1.0;
213:     }
214:     for (ni=1;ni<n && ss[ni]==ss[0]; ni++);
215:     n0 = ni-1;
216:     n1 = n_-ni;
217:     for (j=0;j<n-2;j++) {
218:       PetscBLASIntCast(n-j-1,&m);
219:       /* Forming and applying reflectors */
220:       if (n0 > 1) {
221:         PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n0,A+ni-n0+j*lda,A+ni-n0+j*lda+1,&inc,&tau));
222:         /* Apply reflector */
223:         if (PetscAbsScalar(tau) != 0.0) {
224:           t=*(A+ni-n0+j*lda);  *(A+ni-n0+j*lda)=1.0;
225:           PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&m,&n0,A+ni-n0+j*lda,&inc,&tau,A+j+1+(j+1)*lda,&lda_,work+nwu));
226:           PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n0,&m,A+ni-n0+j*lda,&inc,&tau,A+j+1+(j+1)*lda,&lda_,work+nwu));
227:           /* Update Q */
228:           PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n0,A+ni-n0+j*lda,&inc,&tau,Q+(j+1)*ldq,&ldq_,work+nwu));
229:           *(A+ni-n0+j*lda) = t;
230:           for (i=1;i<n0;i++) {
231:             *(A+ni-n0+j*lda+i) = 0.0;  *(A+j+(ni-n0+i)*lda) = 0.0;
232:           }
233:           *(A+j+(ni-n0)*lda) = *(A+ni-n0+j*lda);
234:         }
235:       }
236:       if (n1 > 1) {
237:         PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n1,A+n-n1+j*lda,A+n-n1+j*lda+1,&inc,&tau));
238:         /* Apply reflector */
239:         if (PetscAbsScalar(tau) != 0.0) {
240:           t=*(A+n-n1+j*lda);  *(A+n-n1+j*lda)=1.0;
241:           PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&m,&n1,A+n-n1+j*lda,&inc,&tau,A+j+1+(n-n1)*lda,&lda_,work+nwu));
242:           PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n1,&m,A+n-n1+j*lda,&inc,&tau,A+n-n1+(j+1)*lda,&lda_,work+nwu));
243:           /* Update Q */
244:           PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n1,A+n-n1+j*lda,&inc,&tau,Q+(n-n1)*ldq,&ldq_,work+nwu));
245:           *(A+n-n1+j*lda) = t;
246:           for (i=1;i<n1;i++) {
247:             *(A+n-n1+i+j*lda) = 0.0;  *(A+j+(n-n1+i)*lda) = 0.0;
248:           }
249:           *(A+j+(n-n1)*lda) = *(A+n-n1+j*lda);
250:         }
251:       }
252:       /* Hyperbolic rotation */
253:       if (n0 > 0 && n1 > 0) {
254:         HRGen(PetscRealPart(A[ni-n0+j*lda]),PetscRealPart(A[n-n1+j*lda]),&type,&cs,&sn,&r,&cond);
255:         /* Check condition number */
256:         if (cond > 1.0/(10*PETSC_SQRT_MACHINE_EPSILON)) {
257:           breakdown = PETSC_TRUE;
258:           k++;
259:           if (k==n || flip)
260:             SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Breakdown in construction of hyperbolic transformation");
261:           break;
262:         }
263:         A[ni-n0+j*lda] = r; A[n-n1+j*lda] = 0.0;
264:         A[j+(ni-n0)*lda] = r; A[j+(n-n1)*lda] = 0.0;
265:         /* Apply to A */
266:         HRApply(m,A+j+1+(ni-n0)*lda,1,A+j+1+(n-n1)*lda,1,cs,-sn);
267:         HRApply(m,A+ni-n0+(j+1)*lda,lda,A+n-n1+(j+1)*lda,lda,cs,-sn);

269:         /* Update Q */
270:         HRApply(n,Q+(ni-n0)*ldq,1,Q+(n-n1)*ldq,1,cs,-sn);
271:         if (type==2) {
272:           ss[ni-n0] = -ss[ni-n0]; ss[n-n1] = -ss[n-n1];
273:           n0++;ni++;n1--;
274:         }
275:       }
276:       if (n0>0) n0--;
277:       else n1--;
278:     }
279:   }

281:   /* flip matrices */
282:   if (flip) {
283:     for (i=0;i<n-1;i++) {
284:       d[i] = PetscRealPart(A[n-i-1+(n-i-1)*lda]);
285:       e[i] = PetscRealPart(A[n-i-1+(n-i-2)*lda]);
286:       s[i] = ss[n-i-1];
287:     }
288:     s[n-1] = ss[0];
289:     d[n-1] = PetscRealPart(A[0]);
290:     for (i=0;i<n;i++) {
291:       ierr=PetscMemcpy(work+i*n,Q+i*ldq,n*sizeof(PetscScalar));
292:     }
293:     for (i=0;i<n;i++)
294:       for (j=0;j<n;j++)
295:         Q[i+j*ldq] = work[i+(n-j-1)*n];
296:   } else {
297:     for (i=0;i<n-1;i++) {
298:       d[i] = PetscRealPart(A[i+i*lda]);
299:       e[i] = PetscRealPart(A[i+1+i*lda]);
300:       s[i] = ss[i];
301:     }
302:     s[n-1] = ss[n-1];
303:     d[n-1] = PetscRealPart(A[n-1 + (n-1)*lda]);
304:   }
305:   return(0);
306: #endif
307: }

311: static PetscErrorCode MadeHRtr(PetscInt sz,PetscInt n,PetscInt idx0,PetscInt n0,PetscInt idx1,PetscInt n1,struct HRtr *tr1,struct HRtr *tr2,PetscReal *ncond,PetscScalar *work)
312: {
313: #if defined(SLEPC_MISSING_LAPACK_LARFG) || defined(SLEPC_MISSING_LAPACK_LARF)
315:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LARFG/LARF - Lapack routines are unavailable");
316: #else
318:   PetscScalar    *x,*y;
319:   PetscReal       ncond2;
320:   PetscBLASInt   n0_,n1_,inc=1;

323:   /* Hyperbolic transformation to make zeros in x */
324:   x = tr1->data;
325:   tr1->n[0] = n0;
326:   tr1->n[1] = n1;
327:   tr1->idx[0] = idx0;
328:   tr1->idx[1] = idx1;
329:   PetscBLASIntCast(tr1->n[0],&n0_);
330:   PetscBLASIntCast(tr1->n[1],&n1_);
331:   if (tr1->n[0] > 1) {
332:     PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n0_,x+tr1->idx[0],x+tr1->idx[0]+1,&inc,tr1->tau));
333:   }
334:   if (tr1->n[1]> 1) {
335:     PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n1_,x+tr1->idx[1],x+tr1->idx[1]+1,&inc,tr1->tau+1));
336:   }
337:   if (tr1->idx[0]<tr1->idx[1]) {
338:     HRGen(PetscRealPart(x[tr1->idx[0]]),PetscRealPart(x[tr1->idx[1]]),&(tr1->type),&(tr1->cs),&(tr1->sn),&(tr1->alpha),ncond);  
339:   } else {
340:     tr1->alpha = PetscRealPart(x[tr1->idx[0]]);
341:     *ncond = 1.0;
342:   }
343:   if (sz==2) {
344:     y = tr2->data;
345:     /* Apply first transformation to second column */
346:     if (tr1->n[0] > 1 && PetscAbsScalar(tr1->tau[0])!=0.0) {
347:       x[tr1->idx[0]] = 1.0;
348:       PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n0_,&inc,x+tr1->idx[0],&inc,tr1->tau,y+tr1->idx[0],&n0_,work));
349:     }
350:     if (tr1->n[1] > 1 && PetscAbsScalar(tr1->tau[1])!=0.0) {
351:       x[tr1->idx[1]] = 1.0;
352:       PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n1_,&inc,x+tr1->idx[1],&inc,tr1->tau+1,y+tr1->idx[1],&n1_,work));
353:     }
354:     if (tr1->idx[0]<tr1->idx[1]) {
355:       HRApply(1,y+tr1->idx[0],1,y+tr1->idx[1],1,tr1->cs,-tr1->sn);
356:     }
357:     tr2->n[0] = tr1->n[0];
358:     tr2->n[1] = tr1->n[1];
359:     tr2->idx[0] = tr1->idx[0];
360:     tr2->idx[1] = tr1->idx[1];
361:     if (tr1->idx[0]<tr1->idx[1] && tr1->type==2) {
362:       tr2->idx[1]++; tr2->n[1]--; tr2->n[0]++;
363:     }
364:     if (tr2->n[0]>0) {
365:       tr2->n[0]--; tr2->idx[0]++;
366:       if (tr2->n[1]==0) tr2->idx[1] = tr2->idx[0];
367:     } else {
368:       tr2->n[1]--; tr2->idx[1]++; tr2->idx[0] = tr2->idx[1];
369:     }
370:     /* Hyperbolic transformation to make zeros in y */
371:     PetscBLASIntCast(tr2->n[0],&n0_);
372:     PetscBLASIntCast(tr2->n[1],&n1_);
373:     if (tr2->n[0] > 1) {
374:       PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n0_,y+tr2->idx[0],y+tr2->idx[0]+1,&inc,tr2->tau));
375:     }
376:     if (tr2->n[1]> 1) {
377:       PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n1_,y+tr2->idx[1],y+tr2->idx[1]+1,&inc,tr2->tau+1));
378:     }
379:     if (tr2->idx[0]<tr2->idx[1]) {
380:       HRGen(PetscRealPart(y[tr2->idx[0]]),PetscRealPart(y[tr2->idx[1]]),&(tr2->type),&(tr2->cs),&(tr2->sn),&(tr2->alpha),&ncond2);  
381:     } else {
382:     tr2->alpha = PetscRealPart(y[tr2->idx[0]]);
383:     ncond2 = 1.0;
384:     }
385:     if (ncond2>*ncond) *ncond = ncond2;
386:   }
387:   return(0);
388: #endif
389: }

393: /*
394:   Auxiliary function to try perform one iteration of hr routine,
395:   checking condition number. If it is < tolD, apply the
396:   transformation to H and R, if not, ok=false and it do nothing
397:   tolE, tolerance to exchange complex pairs to improve conditioning
398: */
399: static PetscErrorCode TryHRIt(PetscInt n,PetscInt j,PetscInt sz,PetscScalar *H,PetscInt ldh,PetscScalar *R,PetscInt ldr,PetscReal *s,PetscBool *exg,PetscBool *ok,PetscInt *n0,PetscInt *n1,PetscInt *idx0,PetscInt *idx1,PetscReal *cond,PetscScalar *work)
400: {
401: #if defined(SLEPC_MISSING_LAPACK_LARF)
403:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LARF - Lapack routine is unavailable");
404: #else
406:   struct HRtr    *tr1,*tr2,tr1_t,tr2_t,tr1_te,tr2_te;
407:   PetscScalar    *x,*y;
408:   PetscReal      ncond,ncond_e;
409:   PetscInt       nwu=0,i,d=1;
410:   PetscBLASInt   n0_,n1_,inc=1,mh,mr,n_,ldr_,ldh_;
411:   PetscReal      tolD = 1e+5;

414:   if (cond) *cond = 1.0;
415:   PetscBLASIntCast(n,&n_);
416:   PetscBLASIntCast(ldr,&ldr_);
417:   PetscBLASIntCast(ldh,&ldh_);
418:   x = work+nwu;
419:   nwu += n;
420:   PetscMemcpy(x,R+j*ldr,n*sizeof(PetscScalar));
421:   *exg = PETSC_FALSE;
422:   *ok = PETSC_TRUE;
423:   tr1_t.data = x;
424:   if (sz==1) {
425:     /* Hyperbolic transformation to make zeros in x */
426:     MadeHRtr(sz,n,*idx0,*n0,*idx1,*n1,&tr1_t,NULL,&ncond,work+nwu);
427:     /* Check condition number to single column*/
428:     if (ncond>tolD) {
429:       *ok = PETSC_FALSE;
430:     }
431:     tr1 = &tr1_t;
432:     tr2 = &tr2_t;    
433:   } else {
434:     y = work+nwu;
435:     nwu += n;
436:     PetscMemcpy(y,R+(j+1)*ldr,n*sizeof(PetscScalar));
437:     tr2_t.data = y;
438:     MadeHRtr(sz,n,*idx0,*n0,*idx1,*n1,&tr1_t,&tr2_t,&ncond,work+nwu);
439:     /* Computing hyperbolic transformations also for exchanged vectors */
440:     tr1_te.data = work+nwu;
441:     nwu += n;
442:     PetscMemcpy(tr1_te.data,R+(j+1)*ldr,n*sizeof(PetscScalar));
443:     tr2_te.data = work+nwu;
444:     nwu += n;
445:     PetscMemcpy(tr2_te.data,R+j*ldr,n*sizeof(PetscScalar));
446:     MadeHRtr(sz,n,*idx0,*n0,*idx1,*n1,&tr1_te,&tr2_te,&ncond_e,work+nwu);
447:     if (ncond > d*ncond_e) {
448:       *exg = PETSC_TRUE;
449:       tr1 = &tr1_te;
450:       tr2 = &tr2_te;
451:       ncond = ncond_e;
452:     } else {
453:       tr1 = &tr1_t;
454:       tr2 = &tr2_t;
455:     }
456:     if (ncond>tolD) *ok = PETSC_FALSE;
457:   }
458:   if (*ok) {
459:     /* Everything is OK, apply transformations to R and H */
460:     /* First column */
461:     if (cond && *cond<ncond) *cond = ncond;
462:     x = tr1->data;
463:     PetscBLASIntCast(tr1->n[0],&n0_);
464:     PetscBLASIntCast(tr1->n[1],&n1_); 
465:     PetscBLASIntCast(n-j-sz,&mr);
466:     if (tr1->n[0] > 1 && PetscAbsScalar(tr1->tau[0])!=0.0) {
467:       x[tr1->idx[0]] = 1.0;
468:       PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n0_,&mr,x+tr1->idx[0],&inc,tr1->tau,R+(j+sz)*ldr+tr1->idx[0],&ldr_,work+nwu));
469:       PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n0_,x+tr1->idx[0],&inc,tr1->tau,H+(tr1->idx[0])*ldh,&ldh_,work+nwu));
470:     }
471:     if (tr1->n[1] > 1 && PetscAbsScalar(tr1->tau[1])!=0.0) {
472:       x[tr1->idx[1]] = 1.0;
473:       PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n1_,&mr,x+tr1->idx[1],&inc,tr1->tau+1,R+(j+sz)*ldr+tr1->idx[1],&ldr_,work+nwu));
474:       PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n1_,x+tr1->idx[1],&inc,tr1->tau+1,H+(tr1->idx[1])*ldh,&ldh_,work+nwu));
475:     }
476:     if (tr1->idx[0]<tr1->idx[1]) {
477:       HRApply(mr,R+(j+sz)*ldr+tr1->idx[0],ldr,R+(j+sz)*ldr+tr1->idx[1],ldr,tr1->cs,-tr1->sn);
478:       if (tr1->type==1) {
479:         HRApply(n,H+(tr1->idx[0])*ldh,1,H+(tr1->idx[1])*ldh,1,tr1->cs,tr1->sn);
480:       } else {
481:         HRApply(n,H+(tr1->idx[0])*ldh,1,H+(tr1->idx[1])*ldh,1,-tr1->cs,-tr1->sn);
482:         s[tr1->idx[0]] = -s[tr1->idx[0]];
483:         s[tr1->idx[1]] = -s[tr1->idx[1]];
484:       }
485:     }
486:     for (i=0;i<tr1->idx[0];i++) *(R+j*ldr+i) = x[i];  
487:     for (i=tr1->idx[0]+1;i<n;i++) *(R+j*ldr+i) = 0.0;
488:     *(R+j*ldr+tr1->idx[0]) = tr1->alpha;
489:     if (sz==2) {
490:       y = tr2->data;
491:       /* Second column */
492:       PetscBLASIntCast(tr2->n[0],&n0_);
493:       PetscBLASIntCast(tr2->n[1],&n1_); 
494:       PetscBLASIntCast(n-j-sz,&mr);
495:       PetscBLASIntCast(n-tr2->idx[0],&mh);
496:       if (tr2->n[0] > 1 && PetscAbsScalar(tr2->tau[0])!=0.0) {
497:         y[tr2->idx[0]] = 1.0;
498:         PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n0_,&mr,y+tr2->idx[0],&inc,tr2->tau,R+(j+2)*ldr+tr2->idx[0],&ldr_,work+nwu));
499:         PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n0_,y+tr2->idx[0],&inc,tr2->tau,H+(tr2->idx[0])*ldh,&ldh_,work+nwu));
500:       }
501:       if (tr2->n[1] > 1 && PetscAbsScalar(tr2->tau[1])!=0.0) {
502:         y[tr2->idx[1]] = 1.0;
503:         PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n1_,&mr,y+tr2->idx[1],&inc,tr2->tau+1,R+(j+2)*ldr+tr2->idx[1],&ldr_,work+nwu));
504:         PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n1_,y+tr2->idx[1],&inc,tr2->tau+1,H+(tr2->idx[1])*ldh,&ldh_,work+nwu));
505:       }
506:       if (tr2->idx[0]<tr2->idx[1]) {
507:         HRApply(mr,R+(j+2)*ldr+tr2->idx[0],ldr,R+(j+2)*ldr+tr2->idx[1],ldr,tr2->cs,-tr2->sn);
508:         if (tr2->type==1) {
509:           HRApply(n,H+(tr2->idx[0])*ldh,1,H+(tr2->idx[1])*ldh,1,tr2->cs,tr2->sn);
510:         } else {
511:           HRApply(n,H+(tr2->idx[0])*ldh,1,H+(tr2->idx[1])*ldh,1,-tr2->cs,-tr2->sn);
512:           s[tr2->idx[0]] = -s[tr2->idx[0]];
513:           s[tr2->idx[1]] = -s[tr2->idx[1]];
514:         }
515:       }
516:       for (i=0;i<tr2->idx[0]-1;i++) *(R+(j+1)*ldr+i) = y[i];
517:       *(R+(j+1)*ldr+tr2->idx[0]-1) = y[tr2->idx[0]-1];
518:       for (i=tr2->idx[0]+1;i<n;i++) *(R+(j+1)*ldr+i) = 0.0;
519:       *(R+(j+1)*ldr+tr2->idx[0]) = tr2->alpha;
520:       *n0 = tr2->n[0];
521:       *n1 = tr2->n[1];
522:       *idx0 = tr2->idx[0];
523:       *idx1 = tr2->idx[1];
524:       if (tr2->idx[0]<tr2->idx[1] && tr2->type==2) {
525:         (*idx1)++; (*n1)--; (*n0)++;
526:       }
527:     } else {
528:       *n0 = tr1->n[0];
529:       *n1 = tr1->n[1];
530:       *idx0 = tr1->idx[0];
531:       *idx1 = tr1->idx[1];
532:       if (tr1->idx[0]<tr1->idx[1] && tr1->type==2) {
533:         (*idx1)++; (*n1)--; (*n0)++;
534:       }      
535:     }
536:     if (*n0>0) {
537:       (*n0)--; (*idx0)++;
538:       if (*n1==0) *idx1 = *idx0;
539:     } else {
540:       (*n1)--; (*idx1)++; *idx0 = *idx1;
541:     }
542:   }
543:   return(0);
544: #endif
545: }

549: /*
550:   compute V = HR whit H s-orthogonal and R upper triangular
551: */
552: static PetscErrorCode PseudoOrthog_HR(PetscInt *nv,PetscScalar *V,PetscInt ldv,PetscReal *s,PetscScalar *R,PetscInt ldr,PetscBLASInt *perm,PetscBLASInt *cmplxEig,PetscBool *breakdown,PetscScalar *work)
553: {
555:   PetscInt       i,j,n,n0,n1,np,idx0,idx1,sz=1,k=0,t1,t2,nwu=0;
556:   PetscScalar    *col1,*col2;
557:   PetscBool      exg=PETSC_FALSE,ok=PETSC_FALSE;

560:   n = *nv;
561:   col1 = work+nwu;
562:   nwu += n;
563:   col2 = work+nwu;
564:   nwu += n;
565:   /* Sort R and s according to sing(s) */
566:   np = 0;
567:   for (i=0;i<n;i++) if (s[i]>0) np++;
568:   if (s[0]>0) n1 = np;
569:   else n1 = n-np;
570:   n0 = 0;
571:   for (i=0;i<n;i++) {
572:     if (s[i]==s[0]) {
573:       s[n0] = s[0];
574:       perm[n0++] = i;
575:     } else perm[n1++] = i;
576:   }
577:   for (i=n0;i<n;i++) s[i] = -s[0];
578:   n1 -= n0;
579:   idx0 = 0;
580:   idx1 = n0;
581:   if (idx1==n) idx1=idx0;
582:   for (i=0;i<n;i++) {
583:     for (j=0;j<n;j++) R[j*ldr+i] = V[j*ldv+perm[i]];
584:   }
585:   /* Initialize H */
586:   for (i=0;i<n;i++) {
587:     PetscMemzero(V+i*ldv,n*sizeof(PetscScalar));
588:     V[perm[i]+i*ldv] = 1.0;
589:   }
590:   for (i=0;i<n;i++) perm[i] = i;
591:   j = 0;
592:   while (j<n-k) {
593:     if (cmplxEig[j]==0) sz=1;
594:     else sz=2;
595:     TryHRIt(n,j,sz,V,ldv,R,ldr,s,&exg,&ok,&n0,&n1,&idx0,&idx1,NULL,work+nwu);
596:     if (ok) {
597:       if (exg) cmplxEig[j] = -cmplxEig[j];
598:       j = j+sz;
599:     } else { /* to be discarded */
600:       k = k+1;
601:       if (cmplxEig[j]==0) {
602:         if (j<n) {
603:           t1 = perm[j];
604:           for (i=j;i<n-1;i++) perm[i] = perm[i+1];
605:           perm[n-1] = t1;
606:           t1 = cmplxEig[j];
607:           for (i=j;i<n-1;i++) cmplxEig[i] = cmplxEig[i+1];
608:           cmplxEig[n-1] = t1;
609:           PetscMemcpy(col1,R+j*ldr,n*sizeof(PetscScalar));
610:           for (i=j;i<n-1;i++) { 
611:             PetscMemcpy(R+i*ldr,R+(i+1)*ldr,n*sizeof(PetscScalar));
612:           }
613:           PetscMemcpy(R+(n-1)*ldr,col1,n*sizeof(PetscScalar));
614:         }
615:       } else {
616:         k = k+1;
617:         if (j<n-1) {
618:           t1 = perm[j];
619:           t2 = perm[j+1];
620:           for (i=j;i<n-2;i++) perm[i] = perm[i+2];
621:           perm[n-2] = t1;
622:           perm[n-1] = t2;
623:           t1 = cmplxEig[j];
624:           t2 = cmplxEig[j+1];
625:           for (i=j;i<n-2;i++) cmplxEig[i] = cmplxEig[i+2];
626:           cmplxEig[n-2] = t1;
627:           cmplxEig[n-1] = t2;
628:           PetscMemcpy(col1,R+j*ldr,n*sizeof(PetscScalar));
629:           PetscMemcpy(col2,R+(j+1)*ldr,n*sizeof(PetscScalar));
630:           for (i=j;i<n-2;i++) {
631:             PetscMemcpy(R+i*ldr,R+(i+2)*ldr,n*sizeof(PetscScalar));
632:           }
633:           PetscMemcpy(R+(n-2)*ldr,col1,n*sizeof(PetscScalar));
634:           PetscMemcpy(R+(n-1)*ldr,col2,n*sizeof(PetscScalar));
635:         }
636:       }
637:     }
638:   }
639:   if (k!=0) {
640:     if (breakdown) *breakdown = PETSC_TRUE;
641:     *nv = n-k;
642:   }
643:   return(0);
644: }

648: PetscErrorCode DSGHIEPOrthogEigenv(DS ds,DSMatType mat,PetscScalar *wr,PetscScalar *wi,PetscBool accum)
649: {
651:   PetscInt       lws,nwus=0,nwui=0,lwi;
652:   PetscInt       off,n,nv,ld,i,ldr,l;
653:   PetscScalar    *W,*X,*R,*ts,zeroS=0.0,oneS=1.0;
654:   PetscReal      *s,vi,vr,tr,*d,*e;
655:   PetscBLASInt   ld_,n_,nv_,*perm,*cmplxEig;

658:   l = ds->l;
659:   n = ds->n-l;
660:   PetscBLASIntCast(n,&n_);
661:   ld = ds->ld;
662:   PetscBLASIntCast(ld,&ld_);
663:   off = l*ld+l;
664:   s = ds->rmat[DS_MAT_D];
665:   if (!ds->compact) {
666:     for (i=l;i<ds->n;i++) s[i] = PetscRealPart(*(ds->mat[DS_MAT_B]+i*ld+i));
667:   }
668:   lws = n*n+7*n;
669:   lwi = 2*n;
670:   DSAllocateWork_Private(ds,lws,0,lwi);
671:   R = ds->work+nwus;
672:   nwus += n*n;
673:   ldr = n;
674:   perm = ds->iwork + nwui;
675:   nwui += n;
676:   cmplxEig = ds->iwork+nwui;
677:   X = ds->mat[mat];
678:   for (i=0;i<n;i++) {
679: #if defined(PETSC_USE_COMPLEX)
680:     vi = PetscImaginaryPart(wr[l+i]);
681: #else
682:     vi = PetscRealPart(wi[l+i]);
683: #endif
684:     if (vi!=0) {
685:       cmplxEig[i] = 1;
686:       cmplxEig[i+1] = 2;
687:       i++;
688:     } else cmplxEig[i] = 0;
689:   }
690:   nv = n;
691:   
692:   /* Perform HR decomposition */
693:   /* Hyperbolic rotators */
694:   PseudoOrthog_HR(&nv,X+off,ld,s+l,R,ldr,perm,cmplxEig,NULL,ds->work+nwus);
695:   /* Sort wr,wi perm */ 
696:   ts = ds->work+nwus;
697:   PetscMemcpy(ts,wr+l,n*sizeof(PetscScalar));
698:   for (i=0;i<n;i++) wr[i+l] = ts[perm[i]];
699: #if !defined(PETSC_USE_COMPLEX)
700:   PetscMemcpy(ts,wi+l,n*sizeof(PetscScalar));
701:   for (i=0;i<n;i++) wi[i+l] = ts[perm[i]];
702: #endif
703:   /* Projected Matrix */
704:   PetscMemzero(ds->rmat[DS_MAT_T]+2*ld,ld*sizeof(PetscReal));
705:   d = ds->rmat[DS_MAT_T];
706:   e = d+ld;
707:   for (i=0;i<nv;i++) {
708:     if (cmplxEig[i]==0) { /* Real */
709:       d[l+i] = PetscRealPart(wr[l+i]*s[l+i]);
710:       e[l+i] = 0.0;
711:     } else {
712:       vr = PetscRealPart(wr[l+i]);
713: #if defined(PETSC_USE_COMPLEX)
714:       vi = PetscImaginaryPart(wr[l+i]);
715: #else
716:       vi = PetscRealPart(wi[l+i]);
717: #endif
718:       if (cmplxEig[i]==-1) vi = -vi;
719:       tr = PetscRealPart((R[i+(i+1)*ldr]/R[i+i*ldr]))*vi;
720:       d[l+i] = (vr-tr)*s[l+i];
721:       d[l+i+1] = (vr+tr)*s[l+i+1];
722:       e[l+i] = PetscRealPart(s[l+i]*(R[(i+1)+(i+1)*ldr]/R[i+i*ldr])*vi);
723:       e[l+i+1] = 0.0;
724:       i++;
725:     }
726:   }
727:   /* accumulate previous Q */
728:   if (accum) {
729:     PetscBLASIntCast(nv,&nv_);
730:     DSAllocateMat_Private(ds,DS_MAT_W);
731:     W = ds->mat[DS_MAT_W];
732:     DSCopyMatrix_Private(ds,DS_MAT_W,DS_MAT_Q);
733:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&nv_,&n_,&oneS,W+off,&ld_,X+off,&ld_,&zeroS,ds->mat[DS_MAT_Q]+off,&ld_));
734:   } else {
735:     PetscMemzero(ds->mat[DS_MAT_Q],ld*ld*sizeof(PetscScalar));
736:     for (i=0;i<ds->l;i++) *(ds->mat[DS_MAT_Q]+i+i*ld) = 1.0;
737:     for (i=0;i<n;i++) { PetscMemcpy(ds->mat[DS_MAT_Q]+off+i*ld,X+off+i*ld,n*sizeof(PetscScalar)); }
738:   }
739:   ds->t = nv+l;
740:   if (!ds->compact) { DSSwitchFormat_GHIEP(ds,PETSC_FALSE); }  
741:   return(0);
742: }

746: /*
747:    Reduce to tridiagonal-diagonal pair by means of TridiagDiag_HHR.
748: */
749: PetscErrorCode DSIntermediate_GHIEP(DS ds)
750: {
752:   PetscInt       i,ld,off;
753:   PetscInt       nwall,nwallr,nwalli;
754:   PetscScalar    *A,*B,*Q;
755:   PetscReal      *d,*e,*s;

758:   ld = ds->ld;
759:   A = ds->mat[DS_MAT_A];
760:   B = ds->mat[DS_MAT_B];
761:   Q = ds->mat[DS_MAT_Q];
762:   d = ds->rmat[DS_MAT_T];
763:   e = ds->rmat[DS_MAT_T]+ld;
764:   s = ds->rmat[DS_MAT_D];
765:   off = ds->l+ds->l*ld;
766:   PetscMemzero(Q,ld*ld*sizeof(PetscScalar));
767:   nwall = ld*ld+ld;
768:   nwallr = ld;
769:   nwalli = ld;
770:   DSAllocateWork_Private(ds,nwall,nwallr,nwalli);
771:   for (i=0;i<ds->n;i++) Q[i+i*ld]=1.0;
772:   for (i=0;i<ds->n-ds->l;i++) *(ds->perm+i)=i;
773:   if (ds->compact) {
774:     if (ds->state < DS_STATE_INTERMEDIATE) {
775:       DSSwitchFormat_GHIEP(ds,PETSC_FALSE);
776:       TridiagDiag_HHR(ds->k-ds->l+1,A+off,ld,s+ds->l,Q+off,ld,PETSC_TRUE,d+ds->l,e+ds->l,ds->perm,ds->work,ds->rwork,ds->iwork);
777:       ds->k = ds->l;
778:       PetscMemzero(d+2*ld+ds->l,(ds->n-ds->l)*sizeof(PetscReal));
779:     }
780:   } else {
781:     if (ds->state < DS_STATE_INTERMEDIATE) {
782:       for (i=0;i<ds->n;i++) s[i] = PetscRealPart(B[i+i*ld]);
783:       TridiagDiag_HHR(ds->n-ds->l,A+off,ld,s+ds->l,Q+off,ld,PETSC_FALSE,d+ds->l,e+ds->l,ds->perm,ds->work,ds->rwork,ds->iwork);
784:       PetscMemzero(d+2*ld,(ds->n)*sizeof(PetscReal));
785:       ds->k = ds->l;
786:       DSSwitchFormat_GHIEP(ds,PETSC_FALSE);
787:     } else {
788:       DSSwitchFormat_GHIEP(ds,PETSC_TRUE);
789:     }
790:   }
791:   return(0);
792: }