Actual source code: dsghiep.c

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

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

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

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

 26: PetscErrorCode DSAllocate_GHIEP(DS ds,PetscInt ld)
 27: {

 31:   DSAllocateMat_Private(ds,DS_MAT_A);
 32:   DSAllocateMat_Private(ds,DS_MAT_B);
 33:   DSAllocateMat_Private(ds,DS_MAT_Q);
 34:   DSAllocateMatReal_Private(ds,DS_MAT_T);
 35:   DSAllocateMatReal_Private(ds,DS_MAT_D);
 36:   PetscFree(ds->perm);
 37:   PetscMalloc1(ld,&ds->perm);
 38:   PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
 39:   return(0);
 40: }

 44: PetscErrorCode DSSwitchFormat_GHIEP(DS ds,PetscBool tocompact)
 45: {
 47:   PetscReal      *T,*S;
 48:   PetscScalar    *A,*B;
 49:   PetscInt       i,n,ld;

 52:   A = ds->mat[DS_MAT_A];
 53:   B = ds->mat[DS_MAT_B];
 54:   T = ds->rmat[DS_MAT_T];
 55:   S = ds->rmat[DS_MAT_D];
 56:   n = ds->n;
 57:   ld = ds->ld;
 58:   if (tocompact) { /* switch from dense (arrow) to compact storage */
 59:     PetscMemzero(T,3*ld*sizeof(PetscReal));
 60:     PetscMemzero(S,ld*sizeof(PetscReal));
 61:     for (i=0;i<n-1;i++) {
 62:       T[i] = PetscRealPart(A[i+i*ld]);
 63:       T[ld+i] = PetscRealPart(A[i+1+i*ld]);
 64:       S[i] = PetscRealPart(B[i+i*ld]);
 65:     }
 66:     T[n-1] = PetscRealPart(A[n-1+(n-1)*ld]);
 67:     S[n-1] = PetscRealPart(B[n-1+(n-1)*ld]);
 68:     for (i=ds->l;i< ds->k;i++) T[2*ld+i] = PetscRealPart(A[ds->k+i*ld]);
 69:   } else { /* switch from compact (arrow) to dense storage */
 70:     PetscMemzero(A,ld*ld*sizeof(PetscScalar));
 71:     PetscMemzero(B,ld*ld*sizeof(PetscScalar));
 72:     for (i=0;i<n-1;i++) {
 73:       A[i+i*ld] = T[i];
 74:       A[i+1+i*ld] = T[ld+i];
 75:       A[i+(i+1)*ld] = T[ld+i];
 76:       B[i+i*ld] = S[i];
 77:     }
 78:     A[n-1+(n-1)*ld] = T[n-1];
 79:     B[n-1+(n-1)*ld] = S[n-1];
 80:     for (i=ds->l;i<ds->k;i++) {
 81:       A[ds->k+i*ld] = T[2*ld+i];
 82:       A[i+ds->k*ld] = T[2*ld+i];
 83:     }
 84:   }
 85:   return(0);
 86: }

 90: PetscErrorCode DSView_GHIEP(DS ds,PetscViewer viewer)
 91: {
 92:   PetscErrorCode    ierr;
 93:   PetscViewerFormat format;
 94:   PetscInt          i,j;
 95:   PetscReal         value;
 96:   const char        *methodname[] = {
 97:                      "HR method",
 98:                      "QR + Inverse Iteration",
 99:                      "QR",
100:                      "DQDS + Inverse Iteration "
101:   };
102:   const int         nmeth=sizeof(methodname)/sizeof(methodname[0]);

105:   PetscViewerGetFormat(viewer,&format);
106:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
107:     if (ds->method>=nmeth) {
108:       PetscViewerASCIIPrintf(viewer,"solving the problem with: INVALID METHOD\n");
109:     } else {
110:       PetscViewerASCIIPrintf(viewer,"solving the problem with: %s\n",methodname[ds->method]);
111:     }
112:     return(0);
113:   }
114:   if (ds->compact) {
115:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
116:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
117:       PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",ds->n,ds->n);
118:       PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",3*ds->n);
119:       PetscViewerASCIIPrintf(viewer,"zzz = [\n");
120:       for (i=0;i<ds->n;i++) {
121:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,i+1,*(ds->rmat[DS_MAT_T]+i));
122:       }
123:       for (i=0;i<ds->n-1;i++) {
124:         if (*(ds->rmat[DS_MAT_T]+ds->ld+i) !=0 && i!=ds->k-1) {
125:           PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+2,i+1,*(ds->rmat[DS_MAT_T]+ds->ld+i));
126:           PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,i+2,*(ds->rmat[DS_MAT_T]+ds->ld+i));
127:         }
128:       }
129:       for (i = ds->l;i<ds->k;i++) {
130:         if (*(ds->rmat[DS_MAT_T]+2*ds->ld+i)) {
131:           PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",ds->k+1,i+1,*(ds->rmat[DS_MAT_T]+2*ds->ld+i));
132:           PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,ds->k+1,*(ds->rmat[DS_MAT_T]+2*ds->ld+i));
133:         }
134:       }
135:       PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_A]);

137:       PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",ds->n,ds->n);
138:       PetscViewerASCIIPrintf(viewer,"omega = zeros(%D,3);\n",3*ds->n);
139:       PetscViewerASCIIPrintf(viewer,"omega = [\n");
140:       for (i=0;i<ds->n;i++) {
141:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,i+1,*(ds->rmat[DS_MAT_D]+i));
142:       }
143:       PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(omega);\n",DSMatName[DS_MAT_B]);

145:     } else {
146:       PetscViewerASCIIPrintf(viewer,"T\n");
147:       for (i=0;i<ds->n;i++) {
148:         for (j=0;j<ds->n;j++) {
149:           if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
150:           else if (i==j+1 || j==i+1) value = *(ds->rmat[DS_MAT_T]+ds->ld+PetscMin(i,j));
151:           else if ((i<ds->k && j==ds->k) || (i==ds->k && j<ds->k)) value = *(ds->rmat[DS_MAT_T]+2*ds->ld+PetscMin(i,j));
152:           else value = 0.0;
153:           PetscViewerASCIIPrintf(viewer," %18.16e ",value);
154:         }
155:         PetscViewerASCIIPrintf(viewer,"\n");
156:       }
157:       PetscViewerASCIIPrintf(viewer,"omega\n");
158:       for (i=0;i<ds->n;i++) {
159:         for (j=0;j<ds->n;j++) {
160:           if (i==j) value = *(ds->rmat[DS_MAT_D]+i);
161:           else value = 0.0;
162:           PetscViewerASCIIPrintf(viewer," %18.16e ",value);
163:         }
164:         PetscViewerASCIIPrintf(viewer,"\n");
165:       }
166:     }
167:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
168:     PetscViewerFlush(viewer);
169:   } else {
170:     DSViewMat(ds,viewer,DS_MAT_A);
171:     DSViewMat(ds,viewer,DS_MAT_B);
172:   }
173:   if (ds->state>DS_STATE_INTERMEDIATE) {
174:     DSViewMat(ds,viewer,DS_MAT_Q);
175:   }
176:   return(0);
177: }

181: static PetscErrorCode DSVectors_GHIEP_Eigen_Some(DS ds,PetscInt *idx,PetscReal *rnorm)
182: {
183: #if defined(SLEPC_MISSING_LAPACK_LAG2)
185:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAG2 - Lapack routine is unavailable");
186: #else
188:   PetscReal      b[4],M[4],d1,d2,s1,s2,e;
189:   PetscReal      scal1,scal2,wr1,wr2,wi,ep,norm;
190:   PetscScalar    *Q,*X,Y[4],alpha,zeroS = 0.0;
191:   PetscInt       k;
192:   PetscBLASInt   two = 2,n_,ld,one=1;
193: #if !defined(PETSC_USE_COMPLEX)
194:   PetscBLASInt   four=4;
195: #endif

198:   X = ds->mat[DS_MAT_X];
199:   Q = ds->mat[DS_MAT_Q];
200:   k = *idx;
201:   PetscBLASIntCast(ds->n,&n_);
202:   PetscBLASIntCast(ds->ld,&ld);
203:   if (k < ds->n-1) {
204:     e = (ds->compact)?*(ds->rmat[DS_MAT_T]+ld+k):PetscRealPart(*(ds->mat[DS_MAT_A]+(k+1)+ld*k));
205:   } else e = 0.0;
206:   if (e == 0.0) {/* Real */
207:     if (ds->state>=DS_STATE_CONDENSED) {
208:       PetscMemcpy(X+k*ld,Q+k*ld,ld*sizeof(PetscScalar));
209:     } else {
210:       PetscMemzero(X+k*ds->ld,ds->ld*sizeof(PetscScalar));
211:       X[k+k*ds->ld] = 1.0;
212:     }
213:     if (rnorm) {
214:       *rnorm = PetscAbsScalar(X[ds->n-1+k*ld]);
215:     }
216:   } else { /* 2x2 block */
217:     if (ds->compact) {
218:       s1 = *(ds->rmat[DS_MAT_D]+k);
219:       d1 = *(ds->rmat[DS_MAT_T]+k);
220:       s2 = *(ds->rmat[DS_MAT_D]+k+1);
221:       d2 = *(ds->rmat[DS_MAT_T]+k+1);
222:     } else {
223:       s1 = PetscRealPart(*(ds->mat[DS_MAT_B]+k*ld+k));
224:       d1 = PetscRealPart(*(ds->mat[DS_MAT_A]+k+k*ld));
225:       s2 = PetscRealPart(*(ds->mat[DS_MAT_B]+(k+1)*ld+k+1));
226:       d2 = PetscRealPart(*(ds->mat[DS_MAT_A]+k+1+(k+1)*ld));
227:     }
228:     M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
229:     b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
230:     ep = LAPACKlamch_("S");
231:     /* Compute eigenvalues of the block */
232:     PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi));
233:     if (wi==0.0)  /* Real eigenvalues */
234:       SETERRQ(PETSC_COMM_SELF,1,"Real block in DSVectors_GHIEP");
235:     else { /* Complex eigenvalues */
236:       if (scal1<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
237:       wr1 /= scal1; wi /= scal1;
238: #if !defined(PETSC_USE_COMPLEX)
239:       if (SlepcAbs(s1*d1-wr1,wi)<SlepcAbs(s2*d2-wr1,wi)) {
240:         Y[0] = wr1-s2*d2; Y[1] = s2*e; Y[2] = wi; Y[3] = 0.0;
241:       } else {
242:         Y[0] = s1*e; Y[1] = wr1-s1*d1; Y[2] = 0.0; Y[3] = wi;
243:       }
244:       norm = BLASnrm2_(&four,Y,&one);
245:       norm = 1/norm;
246:       if (ds->state >= DS_STATE_CONDENSED) {
247:         alpha = norm;
248:         PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&two,&two,&alpha,ds->mat[DS_MAT_Q]+k*ld,&ld,Y,&two,&zeroS,X+k*ld,&ld));
249:         if (rnorm) *rnorm = SlepcAbsEigenvalue(X[ds->n-1+k*ld],X[ds->n-1+(k+1)*ld]);
250:       } else {
251:         PetscMemzero(X+k*ld,2*ld*sizeof(PetscScalar));
252:         X[k*ld+k] = Y[0]*norm; X[k*ld+k+1] = Y[1]*norm;
253:         X[(k+1)*ld+k] = Y[2]*norm; X[(k+1)*ld+k+1] = Y[3]*norm;
254:       }
255: #else
256:       if (SlepcAbs(s1*d1-wr1,wi)<SlepcAbs(s2*d2-wr1,wi)) {
257:         Y[0] = wr1-s2*d2+PETSC_i*wi; Y[1] = s2*e;
258:       } else {
259:         Y[0] = s1*e; Y[1] = wr1-s1*d1+PETSC_i*wi;
260:       }
261:       norm = BLASnrm2_(&two,Y,&one);
262:       norm = 1/norm;
263:       if (ds->state >= DS_STATE_CONDENSED) {
264:         alpha = norm;
265:         PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&two,&alpha,ds->mat[DS_MAT_Q]+k*ld,&ld,Y,&one,&zeroS,X+k*ld,&one));
266:         if (rnorm) *rnorm = PetscAbsScalar(X[ds->n-1+k*ld]);
267:       } else {
268:         PetscMemzero(X+k*ld,2*ld*sizeof(PetscScalar));
269:         X[k*ld+k] = Y[0]*norm; X[k*ld+k+1] = Y[1]*norm;
270:       }
271:       X[(k+1)*ld+k] = PetscConj(X[k*ld+k]); X[(k+1)*ld+k+1] = PetscConj(X[k*ld+k+1]);
272: #endif
273:       (*idx)++;
274:     }
275:   }
276:   return(0);
277: #endif
278: }

282: PetscErrorCode DSVectors_GHIEP(DS ds,DSMatType mat,PetscInt *k,PetscReal *rnorm)
283: {
284:   PetscInt       i;
285:   PetscReal      e;

289:   switch (mat) {
290:     case DS_MAT_X:
291:     case DS_MAT_Y:
292:       if (k) {
293:         DSVectors_GHIEP_Eigen_Some(ds,k,rnorm);
294:       } else {
295:         for (i=0; i<ds->n; i++) {
296:           e = (ds->compact)?*(ds->rmat[DS_MAT_T]+ds->ld+i):PetscRealPart(*(ds->mat[DS_MAT_A]+(i+1)+ds->ld*i));
297:           if (e == 0.0) {/* real */
298:             if (ds->state >= DS_STATE_CONDENSED) {
299:               PetscMemcpy(ds->mat[mat]+i*ds->ld,ds->mat[DS_MAT_Q]+i*ds->ld,ds->ld*sizeof(PetscScalar));
300:             } else {
301:               PetscMemzero(ds->mat[mat]+i*ds->ld,ds->ld*sizeof(PetscScalar));
302:               *(ds->mat[mat]+i+i*ds->ld) = 1.0;
303:             }
304:           } else {
305:             DSVectors_GHIEP_Eigen_Some(ds,&i,rnorm);
306:           }
307:         }
308:       }
309:       break;
310:     case DS_MAT_U:
311:     case DS_MAT_VT:
312:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
313:       break;
314:     default:
315:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
316:   }
317:   return(0);
318: }

322: /*
323:   Extract the eigenvalues contained in the block-diagonal of the indefinite problem.
324:   Only the index range n0..n1 is processed.
325: */
326: PetscErrorCode DSGHIEPComplexEigs(DS ds,PetscInt n0,PetscInt n1,PetscScalar *wr,PetscScalar *wi)
327: {
328: #if defined(SLEPC_MISSING_LAPACK_LAG2)
330:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAG2 - Lapack routine is unavailable");
331: #else
332:   PetscInt     k,ld;
333:   PetscBLASInt two=2;
334:   PetscScalar  *A,*B;
335:   PetscReal    *D,*T;
336:   PetscReal    b[4],M[4],d1,d2,s1,s2,e;
337:   PetscReal    scal1,scal2,ep,wr1,wr2,wi1;

340:   ld = ds->ld;
341:   A = ds->mat[DS_MAT_A];
342:   B = ds->mat[DS_MAT_B];
343:   D = ds->rmat[DS_MAT_D];
344:   T = ds->rmat[DS_MAT_T];
345:   for (k=n0;k<n1;k++) {
346:     if (k < n1-1) {
347:       e = (ds->compact)?T[ld+k]:PetscRealPart(A[(k+1)+ld*k]);
348:     } else {
349:       e = 0.0;
350:     }
351:     if (e==0.0) {
352:       /* real eigenvalue */
353:       wr[k] = (ds->compact)?T[k]/D[k]:A[k+k*ld]/B[k+k*ld];
354: #if !defined(PETSC_USE_COMPLEX)
355:       wi[k] = 0.0 ;
356: #endif
357:     } else {
358:       /* diagonal block */
359:       if (ds->compact) {
360:         s1 = D[k];
361:         d1 = T[k];
362:         s2 = D[k+1];
363:         d2 = T[k+1];
364:       } else {
365:         s1 = PetscRealPart(B[k*ld+k]);
366:         d1 = PetscRealPart(A[k+k*ld]);
367:         s2 = PetscRealPart(B[(k+1)*ld+k+1]);
368:         d2 = PetscRealPart(A[k+1+(k+1)*ld]);
369:       }
370:       M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
371:       b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
372:       ep = LAPACKlamch_("S");
373:       /* Compute eigenvalues of the block */
374:       PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi1));
375:       if (scal1<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
376:       wr[k] = wr1/scal1;
377:       if (wi1==0.0) { /* Real eigenvalues */
378:         if (scal2<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
379:         wr[k+1] = wr2/scal2;
380: #if !defined(PETSC_USE_COMPLEX)
381:         wi[k] = 0.0;
382:         wi[k+1] = 0.0;
383: #endif
384:       } else { /* Complex eigenvalues */
385: #if !defined(PETSC_USE_COMPLEX)
386:         wr[k+1] = wr[k];
387:         wi[k] = wi1/scal1;
388:         wi[k+1] = -wi[k];
389: #else
390:         wr[k] += PETSC_i*wi1/scal1;
391:         wr[k+1] = PetscConj(wr[k]);
392: #endif
393:       }
394:       k++;
395:     }
396:   }
397: #if defined(PETSC_USE_COMPLEX)
398:   if (wi) {
399:     for (k=n0;k<n1;k++) wi[k] = 0.0;
400:   }
401: #endif
402:   return(0);
403: #endif
404: }

408: PetscErrorCode DSSort_GHIEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
409: {
411:   PetscInt       n,i,*perm;
412:   PetscReal      *d,*e,*s;

415: #if !defined(PETSC_USE_COMPLEX)
417: #endif
418:   n = ds->n;
419:   d = ds->rmat[DS_MAT_T];
420:   e = d + ds->ld;
421:   s = ds->rmat[DS_MAT_D];
422:   DSAllocateWork_Private(ds,ds->ld,ds->ld,0);
423:   perm = ds->perm;
424:   if (!rr) {
425:     rr = wr;
426:     ri = wi;
427:   }
428:   DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_TRUE);
429:   if (!ds->compact) { DSSwitchFormat_GHIEP(ds,PETSC_TRUE); }
430:   PetscMemcpy(ds->work,wr,n*sizeof(PetscScalar));
431:   for (i=ds->l;i<n;i++) wr[i] = *(ds->work+perm[i]);
432: #if !defined(PETSC_USE_COMPLEX)
433:   PetscMemcpy(ds->work,wi,n*sizeof(PetscScalar));
434:   for (i=ds->l;i<n;i++) wi[i] = *(ds->work+perm[i]);
435: #endif
436:   PetscMemcpy(ds->rwork,s,n*sizeof(PetscReal));
437:   for (i=ds->l;i<n;i++) s[i] = *(ds->rwork+perm[i]);
438:   PetscMemcpy(ds->rwork,d,n*sizeof(PetscReal));
439:   for (i=ds->l;i<n;i++) d[i] = *(ds->rwork+perm[i]);
440:   PetscMemcpy(ds->rwork,e,(n-1)*sizeof(PetscReal));
441:   PetscMemzero(e+ds->l,(n-1-ds->l)*sizeof(PetscScalar));
442:   for (i=ds->l;i<n-1;i++) {
443:     if (perm[i]<n-1) e[i] = *(ds->rwork+perm[i]);
444:   }
445:   if (!ds->compact) { DSSwitchFormat_GHIEP(ds,PETSC_FALSE); }
446:   DSPermuteColumns_Private(ds,ds->l,n,DS_MAT_Q,perm);
447:   return(0);
448: }


453: /*
454:   Get eigenvectors with inverse iteration.
455:   The system matrix is in Hessenberg form.
456: */
457: PetscErrorCode DSGHIEPInverseIteration(DS ds,PetscScalar *wr,PetscScalar *wi)
458: {
459: #if defined(SLEPC_MISSING_LAPACK_HSEIN)
461:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"HSEIN - Lapack routine is unavailable");
462: #else
464:   PetscInt       i,off;
465:   PetscBLASInt   *select,*infoC,ld,n1,mout,info;
466:   PetscScalar    *A,*B,*H,*X;
467:   PetscReal      *s,*d,*e;
468: #if defined(PETSC_USE_COMPLEX)
469:   PetscInt       j;
470: #endif

473:   PetscBLASIntCast(ds->ld,&ld);
474:   PetscBLASIntCast(ds->n-ds->l,&n1);
475:   DSAllocateWork_Private(ds,ld*ld+2*ld,ld,2*ld);
476:   DSAllocateMat_Private(ds,DS_MAT_W);
477:   A = ds->mat[DS_MAT_A];
478:   B = ds->mat[DS_MAT_B];
479:   H = ds->mat[DS_MAT_W];
480:   s = ds->rmat[DS_MAT_D];
481:   d = ds->rmat[DS_MAT_T];
482:   e = d + ld;
483:   select = ds->iwork;
484:   infoC = ds->iwork + ld;
485:   off = ds->l+ds->l*ld;
486:   if (ds->compact) {
487:     H[off] = d[ds->l]*s[ds->l];
488:     H[off+ld] = e[ds->l]*s[ds->l];
489:     for (i=ds->l+1;i<ds->n-1;i++) {
490:       H[i+(i-1)*ld] = e[i-1]*s[i];
491:       H[i+i*ld] = d[i]*s[i];
492:       H[i+(i+1)*ld] = e[i]*s[i];
493:     }
494:     H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
495:     H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
496:   } else {
497:     s[ds->l] = PetscRealPart(B[off]);
498:     H[off] = A[off]*s[ds->l];
499:     H[off+ld] = A[off+ld]*s[ds->l];
500:     for (i=ds->l+1;i<ds->n-1;i++) {
501:       s[i] = PetscRealPart(B[i+i*ld]);
502:       H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
503:       H[i+i*ld]     = A[i+i*ld]*s[i];
504:       H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
505:     }
506:     s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
507:     H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
508:     H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
509:   }
510:   DSAllocateMat_Private(ds,DS_MAT_X);
511:   X = ds->mat[DS_MAT_X];
512:   for (i=0;i<n1;i++) select[i] = 1;
513: #if !defined(PETSC_USE_COMPLEX)
514:   PetscStackCallBLAS("LAPACKhsein",LAPACKhsein_("R","N","N",select,&n1,H+off,&ld,wr+ds->l,wi+ds->l,NULL,&ld,X+off,&ld,&n1,&mout,ds->work,NULL,infoC,&info));
515: #else
516:   PetscStackCallBLAS("LAPACKhsein",LAPACKhsein_("R","N","N",select,&n1,H+off,&ld,wr+ds->l,NULL,&ld,X+off,&ld,&n1,&mout,ds->work,ds->rwork,NULL,infoC,&info));

518:   /* Separate real and imaginary part of complex eigenvectors */
519:   for (j=ds->l;j<ds->n;j++) {
520:     if (PetscAbsReal(PetscImaginaryPart(wr[j])) > PetscAbsScalar(wr[j])*PETSC_SQRT_MACHINE_EPSILON) {
521:       for (i=ds->l;i<ds->n;i++) {
522:         X[i+(j+1)*ds->ld] = PetscImaginaryPart(X[i+j*ds->ld]);
523:         X[i+j*ds->ld] = PetscRealPart(X[i+j*ds->ld]);
524:       }
525:       j++;
526:     }
527:   }
528: #endif
529:   if (info<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in hsein routine %D",-i);
530:   if (info>0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Convergence error in hsein routine %D",i);
531:   DSGHIEPOrthogEigenv(ds,DS_MAT_X,wr,wi,PETSC_TRUE);
532:   return(0);
533: #endif
534: }


539: /*
540:    Undo 2x2 blocks that have real eigenvalues.
541: */
542: PetscErrorCode DSGHIEPRealBlocks(DS ds)
543: {
544: #if defined(SLEPC_MISSING_LAPACK_LAG2)
546:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAG2 - Lapack routine is unavailable");
547: #else
549:   PetscInt       i;
550:   PetscReal      e,d1,d2,s1,s2,ss1,ss2,t,dd,ss;
551:   PetscReal      maxy,ep,scal1,scal2,snorm;
552:   PetscReal      *T,*D,b[4],M[4],wr1,wr2,wi;
553:   PetscScalar    *A,*B,Y[4],oneS = 1.0,zeroS = 0.0;
554:   PetscBLASInt   m,two=2,ld;
555:   PetscBool      isreal;

558:   PetscBLASIntCast(ds->ld,&ld);
559:   PetscBLASIntCast(ds->n-ds->l,&m);
560:   A = ds->mat[DS_MAT_A];
561:   B = ds->mat[DS_MAT_B];
562:   T = ds->rmat[DS_MAT_T];
563:   D = ds->rmat[DS_MAT_D];
564:   DSAllocateWork_Private(ds,2*m,0,0);
565:   for (i=ds->l;i<ds->n-1;i++) {
566:     e = (ds->compact)?T[ld+i]:PetscRealPart(A[(i+1)+ld*i]);
567:     if (e != 0.0) { /* 2x2 block */
568:       if (ds->compact) {
569:         s1 = D[i];
570:         d1 = T[i];
571:         s2 = D[i+1];
572:         d2 = T[i+1];
573:       } else {
574:         s1 = PetscRealPart(B[i*ld+i]);
575:         d1 = PetscRealPart(A[i*ld+i]);
576:         s2 = PetscRealPart(B[(i+1)*ld+i+1]);
577:         d2 = PetscRealPart(A[(i+1)*ld+i+1]);
578:       }
579:       isreal = PETSC_FALSE;
580:       if (s1==s2) { /* apply a Jacobi rotation to compute the eigendecomposition */
581:         dd = d1-d2;
582:         if (2*PetscAbsReal(e) <= dd) {
583:           t = 2*e/dd;
584:           t = t/(1 + PetscSqrtReal(1+t*t));
585:         } else {
586:           t = dd/(2*e);
587:           ss = (t>=0)?1.0:-1.0;
588:           t = ss/(PetscAbsReal(t)+PetscSqrtReal(1+t*t));
589:         }
590:         Y[0] = 1/PetscSqrtReal(1 + t*t); Y[3] = Y[0]; /* c */
591:         Y[1] = Y[0]*t; Y[2] = -Y[1]; /* s */
592:         wr1 = d1+t*e;
593:         wr2 = d2-t*e;
594:         ss1 = s1; ss2 = s2;
595:         isreal = PETSC_TRUE;
596:       } else {
597:         ss1 = 1.0; ss2 = 1.0,
598:         M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
599:         b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
600:         ep = LAPACKlamch_("S");

602:         /* Compute eigenvalues of the block */
603:         PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi));
604:         if (wi==0.0) { /* Real eigenvalues */
605:           isreal = PETSC_TRUE;
606:           if (scal1<ep||scal2<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
607:           wr1 /= scal1; wr2 /= scal2;
608:           if (PetscAbsReal(s1*d1-wr1)<PetscAbsReal(s2*d2-wr1)) {
609:             Y[0] = wr1-s2*d2;
610:             Y[1] = s2*e;
611:           } else {
612:             Y[0] = s1*e;
613:             Y[1] = wr1-s1*d1;
614:           }
615:           /* normalize with a signature*/
616:           maxy = PetscMax(PetscAbsScalar(Y[0]),PetscAbsScalar(Y[1]));
617:           scal1 = PetscRealPart(Y[0])/maxy; scal2 = PetscRealPart(Y[1])/maxy;
618:           snorm = scal1*scal1*s1 + scal2*scal2*s2;
619:           if (snorm<0) { ss1 = -1.0; snorm = -snorm; }
620:           snorm = maxy*PetscSqrtReal(snorm); Y[0] = Y[0]/snorm; Y[1] = Y[1]/snorm;
621:           if (PetscAbsReal(s1*d1-wr2)<PetscAbsReal(s2*d2-wr2)) {
622:             Y[2] = wr2-s2*d2;
623:             Y[3] = s2*e;
624:           } else {
625:             Y[2] = s1*e;
626:             Y[3] = wr2-s1*d1;
627:           }
628:           maxy = PetscMax(PetscAbsScalar(Y[2]),PetscAbsScalar(Y[3]));
629:           scal1 = PetscRealPart(Y[2])/maxy; scal2 = PetscRealPart(Y[3])/maxy;
630:           snorm = scal1*scal1*s1 + scal2*scal2*s2;
631:           if (snorm<0) { ss2 = -1.0; snorm = -snorm; }
632:           snorm = maxy*PetscSqrtReal(snorm);Y[2] = Y[2]/snorm; Y[3] = Y[3]/snorm;
633:         }
634:         wr1 *= ss1; wr2 *= ss2;
635:       }
636:       if (isreal) {
637:         if (ds->compact) {
638:           D[i] = ss1;
639:           T[i] = wr1;
640:           D[i+1] = ss2;
641:           T[i+1] = wr2;
642:           T[ld+i] = 0.0;
643:         } else {
644:           B[i*ld+i] = ss1;
645:           A[i*ld+i] = wr1;
646:           B[(i+1)*ld+i+1] = ss2;
647:           A[(i+1)*ld+i+1] = wr2;
648:           A[(i+1)+ld*i] = 0.0;
649:           A[i+ld*(i+1)] = 0.0;
650:         }
651:         PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&two,&two,&oneS,ds->mat[DS_MAT_Q]+ds->l+i*ld,&ld,Y,&two,&zeroS,ds->work,&m));
652:         PetscMemcpy(ds->mat[DS_MAT_Q]+ds->l+i*ld,ds->work,m*sizeof(PetscScalar));
653:         PetscMemcpy(ds->mat[DS_MAT_Q]+ds->l+(i+1)*ld,ds->work+m,m*sizeof(PetscScalar));
654:       }
655:       i++;
656:     }
657:   }
658:   return(0);
659: #endif
660: }

664: PetscErrorCode DSSolve_GHIEP_QR_II(DS ds,PetscScalar *wr,PetscScalar *wi)
665: {
666: #if defined(PETSC_MISSING_LAPACK_HSEQR)
668:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"HSEQR - Lapack routine is unavailable");
669: #else
671:   PetscInt       i,off;
672:   PetscBLASInt   n1,ld,one,info,lwork;
673:   PetscScalar    *H,*A,*B,*Q;
674:   PetscReal      *d,*e,*s;
675: #if defined(PETSC_USE_COMPLEX)
676:   PetscInt       j;
677: #endif

680: #if !defined(PETSC_USE_COMPLEX)
682: #endif
683:   one = 1;
684:   PetscBLASIntCast(ds->n-ds->l,&n1);
685:   PetscBLASIntCast(ds->ld,&ld);
686:   off = ds->l + ds->l*ld;
687:   A = ds->mat[DS_MAT_A];
688:   B = ds->mat[DS_MAT_B];
689:   Q = ds->mat[DS_MAT_Q];
690:   d = ds->rmat[DS_MAT_T];
691:   e = ds->rmat[DS_MAT_T] + ld;
692:   s = ds->rmat[DS_MAT_D];
693:   DSAllocateWork_Private(ds,ld*ld,2*ld,ld*2);
694:   lwork = ld*ld;

696:   /* Quick return if possible */
697:   if (n1 == 1) {
698:     *(Q+off) = 1;
699:     if (!ds->compact) {
700:       d[ds->l] = PetscRealPart(A[off]);
701:       s[ds->l] = PetscRealPart(B[off]);
702:     }
703:     wr[ds->l] = d[ds->l]/s[ds->l];
704:     if (wi) wi[ds->l] = 0.0;
705:     return(0);
706:   }
707:   /* Reduce to pseudotriadiagonal form */
708:   DSIntermediate_GHIEP(ds);

710:   /* Compute Eigenvalues (QR)*/
711:   DSAllocateMat_Private(ds,DS_MAT_W);
712:   H = ds->mat[DS_MAT_W];
713:   if (ds->compact) {
714:     H[off] = d[ds->l]*s[ds->l];
715:     H[off+ld] = e[ds->l]*s[ds->l];
716:     for (i=ds->l+1;i<ds->n-1;i++) {
717:       H[i+(i-1)*ld] = e[i-1]*s[i];
718:       H[i+i*ld]     = d[i]*s[i];
719:       H[i+(i+1)*ld] = e[i]*s[i];
720:     }
721:     H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
722:     H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
723:   } else {
724:     s[ds->l] = PetscRealPart(B[off]);
725:     H[off] = A[off]*s[ds->l];
726:     H[off+ld] = A[off+ld]*s[ds->l];
727:     for (i=ds->l+1;i<ds->n-1;i++) {
728:       s[i] = PetscRealPart(B[i+i*ld]);
729:       H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
730:       H[i+i*ld]     = A[i+i*ld]*s[i];
731:       H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
732:     }
733:     s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
734:     H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
735:     H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
736:   }

738: #if !defined(PETSC_USE_COMPLEX)
739:   PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("E","N",&n1,&one,&n1,H+off,&ld,wr+ds->l,wi+ds->l,NULL,&ld,ds->work,&lwork,&info));
740: #else
741:   PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("E","N",&n1,&one,&n1,H+off,&ld,wr+ds->l,NULL,&ld,ds->work,&lwork,&info));

743:   /* Sort to have consecutive conjugate pairs */
744:   for (i=ds->l;i<ds->n;i++) {
745:       j=i+1;
746:       while (j<ds->n && (PetscAbsScalar(wr[i]-PetscConj(wr[j]))>PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON)) j++;
747:       if (j==ds->n) {
748:         if (PetscAbsReal(PetscImaginaryPart(wr[i]))<PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON) wr[i]=PetscRealPart(wr[i]);
749:         else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Found complex without conjugate pair");
750:       } else { /* complex eigenvalue */
751:         wr[j] = wr[i+1];
752:         if (PetscImaginaryPart(wr[i])<0) wr[i] = PetscConj(wr[i]);
753:         wr[i+1] = PetscConj(wr[i]);
754:         i++;
755:       }
756:   }
757: #endif
758:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",&info);
759:   /* Compute Eigenvectors with Inverse Iteration */
760:   DSGHIEPInverseIteration(ds,wr,wi);

762:   /* Recover eigenvalues from diagonal */
763:   DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
764: #if defined(PETSC_USE_COMPLEX)
765:   if (wi) {
766:     for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
767:   }
768: #endif
769:   return(0);
770: #endif
771: }

775: PetscErrorCode DSSolve_GHIEP_QR(DS ds,PetscScalar *wr,PetscScalar *wi)
776: {
777: #if defined(SLEPC_MISSING_LAPACK_GEEVX)
779:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEEVX - Lapack routine is unavailable");
780: #else
782:   PetscInt       i,off,nwu=0,n,lw,lwr,nwru=0;
783:   PetscBLASInt   n_,ld,info,lwork,ilo,ihi;
784:   PetscScalar    *H,*A,*B,*Q,*X;
785:   PetscReal      *d,*s,*scale,nrm,*rcde,*rcdv;
786: #if defined(PETSC_USE_COMPLEX)
787:   PetscInt       j,k;
788: #endif

791: #if !defined(PETSC_USE_COMPLEX)
793: #endif
794:   n = ds->n-ds->l;
795:   PetscBLASIntCast(n,&n_);
796:   PetscBLASIntCast(ds->ld,&ld);
797:   off = ds->l + ds->l*ld;
798:   A = ds->mat[DS_MAT_A];
799:   B = ds->mat[DS_MAT_B];
800:   Q = ds->mat[DS_MAT_Q];
801:   d = ds->rmat[DS_MAT_T];
802:   s = ds->rmat[DS_MAT_D];
803:   lw = 14*ld+ld*ld;
804:   lwr = 7*ld;
805:   DSAllocateWork_Private(ds,lw,lwr,0);
806:   scale = ds->rwork+nwru;
807:   nwru += ld;
808:   rcde = ds->rwork+nwru;
809:   nwru += ld;
810:   rcdv = ds->rwork+nwru;
811:   nwru += ld;
812:   /* Quick return if possible */
813:   if (n_ == 1) {
814:     *(Q+off) = 1;
815:     if (!ds->compact) {
816:       d[ds->l] = PetscRealPart(A[off]);
817:       s[ds->l] = PetscRealPart(B[off]);
818:     }
819:     wr[ds->l] = d[ds->l]/s[ds->l];
820:     if (wi) wi[ds->l] = 0.0;
821:     return(0);
822:   }

824:   /* Form pseudo-symmetric matrix */
825:   H =  ds->work+nwu;
826:   nwu += n*n;
827:   PetscMemzero(H,n*n*sizeof(PetscScalar));
828:   if (ds->compact) {
829:     for (i=0;i<n-1;i++) {
830:       H[i+i*n]     = s[ds->l+i]*d[ds->l+i];
831:       H[i+1+i*n]   = s[ds->l+i+1]*d[ld+ds->l+i];
832:       H[i+(i+1)*n] = s[ds->l+i]*d[ld+ds->l+i];
833:     }
834:     H[n-1+(n-1)*n] = s[ds->l+n-1]*d[ds->l+n-1];
835:     for (i=0;i<ds->k-ds->l;i++) {
836:       H[ds->k-ds->l+i*n] = s[ds->k]*d[2*ld+ds->l+i];
837:       H[i+(ds->k-ds->l)*n] = s[i+ds->l]*d[2*ld+ds->l+i];
838:     }
839:   } else {
840:     for (i=0;i<n-1;i++) {
841:       H[i+i*n]     = B[off+i+i*ld]*A[off+i+i*ld];
842:       H[i+1+i*n]   = B[off+i+1+(i+1)*ld]*A[off+i+1+i*ld];
843:       H[i+(i+1)*n] = B[off+i+i*ld]*A[off+i+(i+1)*ld];
844:     }
845:     H[n-1+(n-1)*n] = B[off+n-1+(n-1)*ld]*A[off+n-1+(n-1)*n];
846:     for (i=0;i<ds->k-ds->l;i++) {
847:       H[ds->k-ds->l+i*n] = B[ds->k*(1+ld)]*A[off+ds->k-ds->l+i*ld];
848:       H[i+(ds->k-ds->l)*n] = B[(i+ds->l)*(1+ld)]*A[off+i+(ds->k-ds->l)*ld];
849:     }
850:   }
851:  
852:   /* Compute eigenpairs */
853:   PetscBLASIntCast(lw-nwu,&lwork);  
854:   DSAllocateMat_Private(ds,DS_MAT_X);
855:   X = ds->mat[DS_MAT_X];
856: #if !defined(PETSC_USE_COMPLEX)
857:   PetscStackCallBLAS("LAPACKgeevx",LAPACKgeevx_("B","N","V","N",&n_,H,&n_,wr+ds->l,wi+ds->l,NULL,&ld,X+off,&ld,&ilo,&ihi,scale,&nrm,rcde,rcdv,ds->work+nwu,&lwork,NULL,&info));
858: #else
859:   PetscStackCallBLAS("LAPACKgeevx",LAPACKgeevx_("B","N","V","N",&n_,H,&n_,wr+ds->l,NULL,&ld,X+off,&ld,&ilo,&ihi,scale,&nrm,rcde,rcdv,ds->work+nwu,&lwork,ds->rwork+nwru,&info));

861:   /* Sort to have consecutive conjugate pairs 
862:      Separate real and imaginary part of complex eigenvectors*/
863:   for (i=ds->l;i<ds->n;i++) {
864:     j=i+1;
865:     while (j<ds->n && (PetscAbsScalar(wr[i]-PetscConj(wr[j]))>PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON)) j++;
866:     if (j==ds->n) {
867:       if (PetscAbsReal(PetscImaginaryPart(wr[i]))<PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON) {
868:         wr[i]=PetscRealPart(wr[i]); /* real eigenvalue */
869:         for (k=ds->l;k<ds->n;k++) {
870:           X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
871:         }
872:       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Found complex without conjugate pair");
873:     } else { /* complex eigenvalue */
874:       if (j!=i+1) {
875:         wr[j] = wr[i+1];
876:         PetscMemcpy(X+j*ds->ld,X+(i+1)*ds->ld,ds->ld*sizeof(PetscScalar));
877:       }
878:       if (PetscImaginaryPart(wr[i])<0) {
879:         wr[i] = PetscConj(wr[i]);
880:         for (k=ds->l;k<ds->n;k++) {
881:           X[k+(i+1)*ds->ld] = -PetscImaginaryPart(X[k+i*ds->ld]);
882:           X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
883:         }
884:       } else {
885:         for (k=ds->l;k<ds->n;k++) {
886:           X[k+(i+1)*ds->ld] = PetscImaginaryPart(X[k+i*ds->ld]);
887:           X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
888:         }
889:       }
890:       wr[i+1] = PetscConj(wr[i]);
891:       i++;
892:     }
893:   }
894: #endif
895:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",&info);

897:   /* Compute real s-orthonormal basis */
898:   DSGHIEPOrthogEigenv(ds,DS_MAT_X,wr,wi,PETSC_FALSE);

900:   /* Recover eigenvalues from diagonal */
901:   DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
902: #if defined(PETSC_USE_COMPLEX)
903:   if (wi) {
904:     for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
905:   }
906: #endif
907:   return(0);
908: #endif
909: }

913: PetscErrorCode DSNormalize_GHIEP(DS ds,DSMatType mat,PetscInt col)
914: {
916:   PetscInt       i,i0,i1;
917:   PetscBLASInt   ld,n,one = 1;
918:   PetscScalar    *A = ds->mat[DS_MAT_A],norm,*x;
919: #if !defined(PETSC_USE_COMPLEX)
920:   PetscScalar    norm0;
921: #endif

924:   switch (mat) {
925:     case DS_MAT_X:
926:     case DS_MAT_Y:
927:     case DS_MAT_Q:
928:       /* Supported matrices */
929:       break;
930:     case DS_MAT_U:
931:     case DS_MAT_VT:
932:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
933:       break;
934:     default:
935:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
936:   }

938:   PetscBLASIntCast(ds->n,&n);
939:   PetscBLASIntCast(ds->ld,&ld);
940:   DSGetArray(ds,mat,&x);
941:   if (col < 0) {
942:     i0 = 0; i1 = ds->n;
943:   } else if (col>0 && A[ds->ld*(col-1)+col] != 0.0) {
944:     i0 = col-1; i1 = col+1;
945:   } else {
946:     i0 = col; i1 = col+1;
947:   }
948:   for (i=i0; i<i1; i++) {
949: #if !defined(PETSC_USE_COMPLEX)
950:     if (i<n-1 && A[ds->ld*i+i+1] != 0.0) {
951:       norm = BLASnrm2_(&n,&x[ld*i],&one);
952:       norm0 = BLASnrm2_(&n,&x[ld*(i+1)],&one);
953:       norm = 1.0/SlepcAbsEigenvalue(norm,norm0);
954:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,&x[ld*i],&one));
955:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,&x[ld*(i+1)],&one));
956:       i++;
957:     } else
958: #endif
959:     {
960:       norm = BLASnrm2_(&n,&x[ld*i],&one);
961:       norm = 1.0/norm;
962:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,&x[ld*i],&one));
963:     }
964:   }
965:   return(0);
966: }

970: PETSC_EXTERN PetscErrorCode DSCreate_GHIEP(DS ds)
971: {
973:   ds->ops->allocate      = DSAllocate_GHIEP;
974:   ds->ops->view          = DSView_GHIEP;
975:   ds->ops->vectors       = DSVectors_GHIEP;
976:   ds->ops->solve[0]      = DSSolve_GHIEP_HZ;
977:   ds->ops->solve[1]      = DSSolve_GHIEP_QR_II;
978:   ds->ops->solve[2]      = DSSolve_GHIEP_QR;
979:   ds->ops->solve[3]      = DSSolve_GHIEP_DQDS_II;
980:   ds->ops->sort          = DSSort_GHIEP;
981:   ds->ops->normalize     = DSNormalize_GHIEP;
982:   return(0);
983: }