Actual source code: dsnhep.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:       
  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: */

 22: #include <slepc-private/dsimpl.h>      /*I "slepcds.h" I*/
 23: #include <slepcblaslapack.h>

 27: PetscErrorCode DSAllocate_NHEP(DS ds,PetscInt ld)
 28: {

 32:   DSAllocateMat_Private(ds,DS_MAT_A);
 33:   DSAllocateMat_Private(ds,DS_MAT_Q);
 34:   PetscFree(ds->perm);
 35:   PetscMalloc(ld*sizeof(PetscInt),&ds->perm);
 36:   PetscLogObjectMemory(ds,ld*sizeof(PetscInt));
 37:   return(0);
 38: }

 42: PetscErrorCode DSView_NHEP(DS ds,PetscViewer viewer)
 43: {

 47:   DSViewMat_Private(ds,viewer,DS_MAT_A);
 48:   if (ds->state>DS_STATE_INTERMEDIATE) {
 49:     DSViewMat_Private(ds,viewer,DS_MAT_Q);
 50:   }
 51:   if (ds->mat[DS_MAT_X]) {
 52:     DSViewMat_Private(ds,viewer,DS_MAT_X);
 53:   }
 54:   if (ds->mat[DS_MAT_Y]) {
 55:     DSViewMat_Private(ds,viewer,DS_MAT_Y);
 56:   }
 57:   return(0);
 58: }

 62: PetscErrorCode DSVectors_NHEP_Refined_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
 63: {
 64: #if defined(SLEPC_MISSING_LAPACK_GESVD)
 66:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable");
 67: #else
 69:   PetscInt       i,j;
 70:   PetscBLASInt   info,ld,n,n1,lwork,inc=1;
 71:   PetscScalar    sdummy,done=1.0,zero=0.0;
 72:   PetscReal      *sigma;
 73:   PetscBool      iscomplex = PETSC_FALSE;
 74:   PetscScalar    *A = ds->mat[DS_MAT_A];
 75:   PetscScalar    *Q = ds->mat[DS_MAT_Q];
 76:   PetscScalar    *X = ds->mat[left?DS_MAT_Y:DS_MAT_X];
 77:   PetscScalar    *W;

 80:   if (left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for left vectors");
 81:   n  = PetscBLASIntCast(ds->n);
 82:   ld = PetscBLASIntCast(ds->ld);
 83:   n1 = n+1;
 84:   if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
 85:   if (iscomplex) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for complex eigenvalues yet");
 86:   DSAllocateWork_Private(ds,5*ld,6*ld,0);
 87:   DSAllocateMat_Private(ds,DS_MAT_W);
 88:   W = ds->mat[DS_MAT_W];
 89:   lwork = 5*ld;
 90:   sigma = ds->rwork+5*ld;

 92:   /* build A-w*I in W */
 93:   for (j=0;j<n;j++)
 94:     for (i=0;i<=n;i++)
 95:       W[i+j*ld] = A[i+j*ld];
 96:   for (i=0;i<n;i++)
 97:     W[i+i*ld] -= A[(*k)+(*k)*ld];

 99:   /* compute SVD of W */
100: #if !defined(PETSC_USE_COMPLEX)
101:   LAPACKgesvd_("N","O",&n1,&n,W,&ld,sigma,&sdummy,&ld,&sdummy,&ld,ds->work,&lwork,&info);
102: #else
103:   LAPACKgesvd_("N","O",&n1,&n,W,&ld,sigma,&sdummy,&ld,&sdummy,&ld,ds->work,&lwork,ds->rwork,&info);
104: #endif
105:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESVD %d",info);

107:   /* the smallest singular value is the new error estimate */
108:   if (rnorm) *rnorm = sigma[n-1];

110:   /* update vector with right singular vector associated to smallest singular value,
111:      accumulating the transformation matrix Q */
112:   BLASgemv_("N",&n,&n,&done,Q,&ld,W+n-1,&ld,&zero,X+(*k)*ld,&inc);
113:   return(0);
114: #endif
115: }

119: PetscErrorCode DSVectors_NHEP_Refined_All(DS ds,PetscBool left)
120: {
122:   PetscInt       i;

125:   for (i=0;i<ds->n;i++) {
126:     DSVectors_NHEP_Refined_Some(ds,&i,PETSC_NULL,left);
127:   }
128:   return(0);
129: }

133: PetscErrorCode DSVectors_NHEP_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
134: {
135: #if defined(SLEPC_MISSING_LAPACK_TREVC)
137:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable");
138: #else
140:   PetscInt       i;
141:   PetscBLASInt   mm=1,mout,info,ld,n,inc = 1;
142:   PetscScalar    tmp,done=1.0,zero=0.0;
143:   PetscReal      norm;
144:   PetscBool      iscomplex = PETSC_FALSE;
145:   PetscBLASInt   *select;
146:   PetscScalar    *A = ds->mat[DS_MAT_A];
147:   PetscScalar    *Q = ds->mat[DS_MAT_Q];
148:   PetscScalar    *X = ds->mat[left?DS_MAT_Y:DS_MAT_X];
149:   PetscScalar    *Y;

152:   n  = PetscBLASIntCast(ds->n);
153:   ld = PetscBLASIntCast(ds->ld);
154:   DSAllocateWork_Private(ds,0,0,ld);
155:   select = ds->iwork;
156:   for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;

158:   /* Compute k-th eigenvector Y of A */
159:   Y = X+(*k)*ld;
160:   select[*k] = (PetscBLASInt)PETSC_TRUE;
161: #if !defined(PETSC_USE_COMPLEX)
162:   if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
163:   mm = iscomplex? 2: 1;
164:   if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
165:   DSAllocateWork_Private(ds,3*ld,0,0);
166:   LAPACKtrevc_(left?"L":"R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,&info);
167: #else
168:   DSAllocateWork_Private(ds,2*ld,ld,0);
169:   LAPACKtrevc_(left?"L":"R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,ds->rwork,&info);
170: #endif
171:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xTREVC %d",info);
172:   if (mout != mm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Inconsistent arguments");

174:   /* accumulate and normalize eigenvectors */
175:   if (ds->state>=DS_STATE_CONDENSED) {
176:     PetscMemcpy(ds->work,Y,mout*ld*sizeof(PetscScalar));
177:     BLASgemv_("N",&n,&n,&done,Q,&ld,ds->work,&inc,&zero,Y,&inc);
178: #if !defined(PETSC_USE_COMPLEX)
179:     if (iscomplex) BLASgemv_("N",&n,&n,&done,Q,&ld,ds->work+ld,&inc,&zero,Y+ld,&inc);
180: #endif
181:     norm = BLASnrm2_(&n,Y,&inc);
182: #if !defined(PETSC_USE_COMPLEX)
183:     if (iscomplex) {
184:       tmp = BLASnrm2_(&n,Y+ld,&inc);
185:       norm = SlepcAbsEigenvalue(norm,tmp);
186:     }
187: #endif
188:     tmp = 1.0 / norm;
189:     BLASscal_(&n,&tmp,Y,&inc);
190: #if !defined(PETSC_USE_COMPLEX)
191:     if (iscomplex) BLASscal_(&n,&tmp,Y+ld,&inc);
192: #endif
193:   }

195:   /* set output arguments */
196:   if (iscomplex) (*k)++;
197:   if (rnorm) {
198:     if (iscomplex) *rnorm = SlepcAbsEigenvalue(Y[n-1],Y[n-1+ld]);
199:     else *rnorm = PetscAbsScalar(Y[n-1]);
200:   }
201:   return(0);
202: #endif
203: }

207: PetscErrorCode DSVectors_NHEP_Eigen_All(DS ds,PetscBool left)
208: {
209: #if defined(SLEPC_MISSING_LAPACK_TREVC)
211:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable");
212: #else
214:   PetscBLASInt   n,ld,mout,info;
215:   PetscScalar    *X,*Y,*A = ds->mat[DS_MAT_A];
216:   const char     *side,*back;

219:   n  = PetscBLASIntCast(ds->n);
220:   ld = PetscBLASIntCast(ds->ld);
221:   if (left) {
222:     X = PETSC_NULL;
223:     Y = ds->mat[DS_MAT_Y];
224:     side = "L";
225:   } else {
226:     X = ds->mat[DS_MAT_X];
227:     Y = PETSC_NULL;
228:     side = "R";
229:   }
230:   if (ds->state>=DS_STATE_CONDENSED) {
231:     /* DSSolve() has been called, backtransform with matrix Q */
232:     back = "B";
233:     PetscMemcpy(left?Y:X,ds->mat[DS_MAT_Q],ld*ld*sizeof(PetscScalar));
234:   } else back = "A";
235: #if !defined(PETSC_USE_COMPLEX)
236:   DSAllocateWork_Private(ds,3*ld,0,0);
237:   LAPACKtrevc_(side,back,PETSC_NULL,&n,A,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,&info);
238: #else
239:   DSAllocateWork_Private(ds,2*ld,ld,0);
240:   LAPACKtrevc_(side,back,PETSC_NULL,&n,A,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info);
241: #endif
242:   if (info) SETERRQ1(((PetscObject)ds)->comm,PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);
243:   return(0);
244: #endif
245: }

249: PetscErrorCode DSVectors_NHEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
250: {

254:   switch (mat) {
255:     case DS_MAT_X:
256:       if (ds->refined) {
257:         if (!ds->extrarow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Refined vectors require activating the extra row");
258:         if (j) {
259:           DSVectors_NHEP_Refined_Some(ds,j,rnorm,PETSC_FALSE);
260:         } else {
261:           DSVectors_NHEP_Refined_All(ds,PETSC_FALSE);
262:         }
263:       } else {
264:         if (j) {
265:           DSVectors_NHEP_Eigen_Some(ds,j,rnorm,PETSC_FALSE);
266:         } else {
267:           DSVectors_NHEP_Eigen_All(ds,PETSC_FALSE);
268:         }
269:       }
270:       break;
271:     case DS_MAT_Y:
272:       if (ds->refined) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
273:       if (j) {
274:         DSVectors_NHEP_Eigen_Some(ds,j,rnorm,PETSC_TRUE);
275:       } else {
276:         DSVectors_NHEP_Eigen_All(ds,PETSC_TRUE);
277:       }
278:       break;
279:     case DS_MAT_U:
280:     case DS_MAT_VT:
281:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
282:       break;
283:     default:
284:       SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
285:   }
286:   return(0);
287: }

291: PetscErrorCode DSNormalize_NHEP(DS ds,DSMatType mat,PetscInt col)
292: {
294:   PetscInt       i,i0,i1;
295:   PetscBLASInt   ld,n,one = 1;
296:   PetscScalar    *A = ds->mat[DS_MAT_A],norm,*x;
297: #if !defined(PETSC_USE_COMPLEX)
298:   PetscScalar    norm0;
299: #endif

302:   switch (mat) {
303:     case DS_MAT_X:
304:     case DS_MAT_Y:
305:     case DS_MAT_Q:
306:       /* Supported matrices */
307:       break;
308:     case DS_MAT_U:
309:     case DS_MAT_VT:
310:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
311:       break;
312:     default:
313:       SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
314:   }

316:   n  = PetscBLASIntCast(ds->n);
317:   ld = PetscBLASIntCast(ds->ld);
318:   DSGetArray(ds,mat,&x);
319:   if (col < 0) {
320:     i0 = 0; i1 = ds->n;
321:   } else if(col>0 && A[ds->ld*(col-1)+col] != 0.0) {
322:     i0 = col-1; i1 = col+1;
323:   } else {
324:     i0 = col; i1 = col+1;
325:   }
326:   for(i=i0; i<i1; i++) {
327: #if !defined(PETSC_USE_COMPLEX)
328:     if(i<n-1 && A[ds->ld*i+i+1] != 0.0) {
329:       norm = BLASnrm2_(&n,&x[ld*i],&one);
330:       norm0 = BLASnrm2_(&n,&x[ld*(i+1)],&one);
331:       norm = 1.0/SlepcAbsEigenvalue(norm,norm0);
332:       BLASscal_(&n,&norm,&x[ld*i],&one);
333:       BLASscal_(&n,&norm,&x[ld*(i+1)],&one);
334:       i++;
335:     } else
336: #endif
337:     {
338:       norm = BLASnrm2_(&n,&x[ld*i],&one);
339:       norm = 1.0/norm;
340:       BLASscal_(&n,&norm,&x[ld*i],&one);
341:      }
342:   }
343:   return(0);
344: }

348: PetscErrorCode DSSort_NHEP_Arbitrary(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
349: {
350: #if defined(SLEPC_MISSING_LAPACK_TRSEN)
352:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TRSEN - Lapack routine is unavailable");
353: #else
355:   PetscInt       i;
356:   PetscBLASInt   info,n,ld,mout,lwork,*selection;
357:   PetscScalar    *T = ds->mat[DS_MAT_A],*Q = ds->mat[DS_MAT_Q],*work;
358: #if !defined(PETSC_USE_COMPLEX)
359:   PetscBLASInt   *iwork,liwork;
360: #endif

363:   if (!k) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must supply argument k");
364:   n  = PetscBLASIntCast(ds->n);
365:   ld = PetscBLASIntCast(ds->ld);
366: #if !defined(PETSC_USE_COMPLEX)
367:   lwork = n;
368:   liwork = 1;
369:   DSAllocateWork_Private(ds,lwork,0,liwork+n);
370:   work = ds->work;
371:   lwork = ds->lwork;
372:   selection = ds->iwork;
373:   iwork = ds->iwork + n;
374:   liwork = ds->liwork - n;
375: #else
376:   lwork = 1;
377:   DSAllocateWork_Private(ds,lwork,0,n);
378:   work = ds->work;
379:   selection = ds->iwork;
380: #endif
381:   /* Compute the selected eigenvalue to be in the leading position */
382:   DSSortEigenvalues_Private(ds,rr,ri,ds->perm,PETSC_FALSE);
383:   PetscMemzero(selection,n*sizeof(PetscBLASInt));
384:   for (i=0;i<*k;i++) selection[ds->perm[i]] = 1;
385: #if !defined(PETSC_USE_COMPLEX)
386:   LAPACKtrsen_("N","V",selection,&n,T,&ld,Q,&ld,wr,wi,&mout,PETSC_NULL,PETSC_NULL,work,&lwork,iwork,&liwork,&info);
387: #else
388:   LAPACKtrsen_("N","V",selection,&n,T,&ld,Q,&ld,wr,&mout,PETSC_NULL,PETSC_NULL,work,&lwork,&info);
389: #endif
390:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xTRSEN %d",info);
391:   *k = mout;
392:   return(0);
393: #endif 
394: }

398: PetscErrorCode DSSort_NHEP_Total(DS ds,PetscScalar *wr,PetscScalar *wi)
399: {
400: #if defined(SLEPC_MISSING_LAPACK_TREXC)
402:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable");
403: #else
405:   PetscScalar    re;
406:   PetscInt       i,j,pos,result;
407:   PetscBLASInt   ifst,ilst,info,n,ld;
408:   PetscScalar    *T = ds->mat[DS_MAT_A];
409:   PetscScalar    *Q = ds->mat[DS_MAT_Q];
410: #if !defined(PETSC_USE_COMPLEX)
411:   PetscScalar    *work,im;
412: #endif

415:   n  = PetscBLASIntCast(ds->n);
416:   ld = PetscBLASIntCast(ds->ld);
417: #if !defined(PETSC_USE_COMPLEX)
418:   DSAllocateWork_Private(ds,ld,0,0);
419:   work = ds->work;
420: #endif
421:   /* selection sort */
422:   for (i=ds->l;i<n-1;i++) {
423:     re = wr[i];
424: #if !defined(PETSC_USE_COMPLEX)
425:     im = wi[i];
426: #endif
427:     pos = 0;
428:     j=i+1; /* j points to the next eigenvalue */
429: #if !defined(PETSC_USE_COMPLEX)
430:     if (im != 0) j=i+2;
431: #endif
432:     /* find minimum eigenvalue */
433:     for (;j<n;j++) {
434: #if !defined(PETSC_USE_COMPLEX)
435:       (*ds->comp_fun)(re,im,wr[j],wi[j],&result,ds->comp_ctx);
436: #else
437:       (*ds->comp_fun)(re,0.0,wr[j],0.0,&result,ds->comp_ctx);
438: #endif
439:       if (result > 0) {
440:         re = wr[j];
441: #if !defined(PETSC_USE_COMPLEX)
442:         im = wi[j];
443: #endif
444:         pos = j;
445:       }
446: #if !defined(PETSC_USE_COMPLEX)
447:       if (wi[j] != 0) j++;
448: #endif
449:     }
450:     if (pos) {
451:       /* interchange blocks */
452:       ifst = PetscBLASIntCast(pos+1);
453:       ilst = PetscBLASIntCast(i+1);
454: #if !defined(PETSC_USE_COMPLEX)
455:       LAPACKtrexc_("V",&n,T,&ld,Q,&ld,&ifst,&ilst,work,&info);
456: #else
457:       LAPACKtrexc_("V",&n,T,&ld,Q,&ld,&ifst,&ilst,&info);
458: #endif
459:       if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
460:       /* recover original eigenvalues from T matrix */
461:       for (j=i;j<n;j++) {
462:         wr[j] = T[j+j*ld];
463: #if !defined(PETSC_USE_COMPLEX)
464:         if (j<n-1 && T[j+1+j*ld] != 0.0) {
465:           /* complex conjugate eigenvalue */
466:           wi[j] = PetscSqrtReal(PetscAbsReal(T[j+1+j*ld])) *
467:                   PetscSqrtReal(PetscAbsReal(T[j+(j+1)*ld]));
468:           wr[j+1] = wr[j];
469:           wi[j+1] = -wi[j];
470:           j++;
471:         } else {
472:           wi[j] = 0.0;
473:         }
474: #endif
475:       }
476:     }
477: #if !defined(PETSC_USE_COMPLEX)
478:     if (wi[i] != 0) i++;
479: #endif
480:   }
481:   return(0);
482: #endif 
483: }

487: PetscErrorCode DSSort_NHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
488: {

492:   if (!rr || wr == rr) {
493:     DSSort_NHEP_Total(ds,wr,wi);
494:   } else {
495:     DSSort_NHEP_Arbitrary(ds,wr,wi,rr,ri,k);
496:   }
497:   return(0);
498: }

502: PetscErrorCode DSUpdateExtraRow_NHEP(DS ds)
503: {
505:   PetscInt       i;
506:   PetscBLASInt   n,ld,incx=1;
507:   PetscScalar    *A,*Q,*x,*y,one=1.0,zero=0.0;

510:   n  = PetscBLASIntCast(ds->n);
511:   ld = PetscBLASIntCast(ds->ld);
512:   A  = ds->mat[DS_MAT_A];
513:   Q  = ds->mat[DS_MAT_Q];
514:   DSAllocateWork_Private(ds,2*ld,0,0);
515:   x = ds->work;
516:   y = ds->work+ld;
517:   for (i=0;i<n;i++) x[i] = A[n+i*ld];
518:   BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx);
519:   for (i=0;i<n;i++) A[n+i*ld] = y[i];
520:   ds->k = n;
521:   return(0);
522: }

526: PetscErrorCode DSSolve_NHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
527: {
528: #if defined(SLEPC_MISSING_LAPACK_GEHRD) || defined(SLEPC_MISSING_LAPACK_ORGHR) || defined(PETSC_MISSING_LAPACK_HSEQR)
530:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEHRD/ORGHR/HSEQR - Lapack routines are unavailable");
531: #else
533:   PetscScalar    *work,*tau;
534:   PetscInt       i,j;
535:   PetscBLASInt   ilo,lwork,info,n,ld;
536:   PetscScalar    *A = ds->mat[DS_MAT_A];
537:   PetscScalar    *Q = ds->mat[DS_MAT_Q];

540: #if !defined(PETSC_USE_COMPLEX)
542: #endif
543:   n   = PetscBLASIntCast(ds->n);
544:   ld  = PetscBLASIntCast(ds->ld);
545:   ilo = PetscBLASIntCast(ds->l+1);
546:   DSAllocateWork_Private(ds,ld+ld*ld,0,0);
547:   tau  = ds->work;
548:   work = ds->work+ld;
549:   lwork = ld*ld;

551:   /* initialize orthogonal matrix */
552:   PetscMemzero(Q,ld*ld*sizeof(PetscScalar));
553:   for (i=0;i<n;i++)
554:     Q[i+i*ld] = 1.0;
555:   if (n==1) return(0);

557:   /* reduce to upper Hessenberg form */
558:   if (ds->state<DS_STATE_INTERMEDIATE) {
559:     LAPACKgehrd_(&n,&ilo,&n,A,&ld,tau,work,&lwork,&info);
560:     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEHRD %d",info);
561:     for (j=0;j<n-1;j++) {
562:       for (i=j+2;i<n;i++) {
563:         Q[i+j*ld] = A[i+j*ld];
564:         A[i+j*ld] = 0.0;
565:       }
566:     }
567:     LAPACKorghr_(&n,&ilo,&n,Q,&ld,tau,work,&lwork,&info);
568:     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGHR %d",info);
569:   }

571:   /* compute the (real) Schur form */
572: #if !defined(PETSC_USE_COMPLEX)
573:   LAPACKhseqr_("S","V",&n,&ilo,&n,A,&ld,wr,wi,Q,&ld,work,&lwork,&info);
574:   for (j=0;j<ds->l;j++) {
575:     if (j==n-1 || A[j+1+j*ld] == 0.0) {
576:       /* real eigenvalue */
577:       wr[j] = A[j+j*ld];
578:       wi[j] = 0.0;
579:     } else {
580:       /* complex eigenvalue */
581:       wr[j] = A[j+j*ld];
582:       wr[j+1] = A[j+j*ld];
583:       wi[j] = PetscSqrtReal(PetscAbsReal(A[j+1+j*ld])) *
584:               PetscSqrtReal(PetscAbsReal(A[j+(j+1)*ld]));
585:       wi[j+1] = -wi[j];
586:       j++;
587:     }
588:   }
589: #else
590:   LAPACKhseqr_("S","V",&n,&ilo,&n,A,&ld,wr,Q,&ld,work,&lwork,&info);
591:   if (wi) for (i=ds->l;i<n;i++) wi[i] = 0.0;
592: #endif
593:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",info);
594:   return(0);
595: #endif
596: }

600: PetscErrorCode DSTruncate_NHEP(DS ds,PetscInt n)
601: {
602:   PetscInt       i,newn,ld=ds->ld,l=ds->l;
603:   PetscScalar    *A;

606:   if (ds->state==DS_STATE_CONDENSED) ds->t = ds->n;
607:   A = ds->mat[DS_MAT_A];
608:   /* be careful not to break a diagonal 2x2 block */
609:   if (A[n+(n-1)*ld]==0.0) newn = n;
610:   else {
611:     if (n<ds->n-1) newn = n+1;
612:     else newn = n-1;
613:   }
614:   if (ds->extrarow && ds->k==ds->n) {
615:     /* copy entries of extra row to the new position, then clean last row */
616:     for (i=l;i<newn;i++) A[newn+i*ld] = A[ds->n+i*ld];
617:     for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
618:   }
619:   ds->k = 0;
620:   ds->n = newn;
621:   return(0);
622: }

626: PetscErrorCode DSCond_NHEP(DS ds,PetscReal *cond)
627: {
628: #if defined(PETSC_MISSING_LAPACK_GETRF) || defined(SLEPC_MISSING_LAPACK_GETRI) || defined(SLEPC_MISSING_LAPACK_LANGE) || defined(SLEPC_MISSING_LAPACK_LANHS)
630:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF/GETRI/LANGE/LANHS - Lapack routines are unavailable");
631: #else
633:   PetscScalar    *work;
634:   PetscReal      *rwork;
635:   PetscBLASInt   *ipiv;
636:   PetscBLASInt   lwork,info,n,ld;
637:   PetscReal      hn,hin;
638:   PetscScalar    *A;

641:   n  = PetscBLASIntCast(ds->n);
642:   ld = PetscBLASIntCast(ds->ld);
643:   lwork = 8*ld;
644:   DSAllocateWork_Private(ds,lwork,ld,ld);
645:   work  = ds->work;
646:   rwork = ds->rwork;
647:   ipiv  = ds->iwork;

649:   /* use workspace matrix W to avoid overwriting A */
650:   DSAllocateMat_Private(ds,DS_MAT_W);
651:   A = ds->mat[DS_MAT_W];
652:   PetscMemcpy(A,ds->mat[DS_MAT_A],sizeof(PetscScalar)*ds->ld*ds->ld);

654:   /* norm of A */
655:   if (ds->state<DS_STATE_INTERMEDIATE) hn = LAPACKlange_("I",&n,&n,A,&ld,rwork);
656:   else hn = LAPACKlanhs_("I",&n,A,&ld,rwork);

658:   /* norm of inv(A) */
659:   LAPACKgetrf_(&n,&n,A,&ld,ipiv,&info);
660:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGETRF %d",info);
661:   LAPACKgetri_(&n,A,&ld,ipiv,work,&lwork,&info);
662:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGETRI %d",info);
663:   hin = LAPACKlange_("I",&n,&n,A,&ld,rwork);

665:   *cond = hn*hin;
666:   return(0);
667: #endif
668: }

672: PetscErrorCode DSTranslateHarmonic_NHEP(DS ds,PetscScalar tau,PetscReal beta,PetscBool recover,PetscScalar *gin,PetscReal *gamma)
673: {
674: #if defined(PETSC_MISSING_LAPACK_GETRF) || defined(PETSC_MISSING_LAPACK_GETRS) 
676:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF/GETRS - Lapack routines are unavailable");
677: #else
679:   PetscInt       i,j;
680:   PetscBLASInt   *ipiv,info,n,ld,one=1,ncol;
681:   PetscScalar    *A,*B,*Q,*g=gin,*ghat;
682:   PetscScalar    done=1.0,dmone=-1.0,dzero=0.0;
683:   PetscReal      gnorm;

686:   n  = PetscBLASIntCast(ds->n);
687:   ld = PetscBLASIntCast(ds->ld);
688:   A  = ds->mat[DS_MAT_A];

690:   if (!recover) {

692:     DSAllocateWork_Private(ds,0,0,ld);
693:     ipiv = ds->iwork;
694:     if (!g) {
695:       DSAllocateWork_Private(ds,ld,0,0);
696:       g = ds->work;
697:     }
698:     /* use workspace matrix W to factor A-tau*eye(n) */
699:     DSAllocateMat_Private(ds,DS_MAT_W);
700:     B = ds->mat[DS_MAT_W];
701:     PetscMemcpy(B,A,sizeof(PetscScalar)*ld*ld);

703:     /* Vector g initialy stores b = beta*e_n^T */
704:     PetscMemzero(g,n*sizeof(PetscScalar));
705:     g[n-1] = beta;
706: 
707:     /* g = (A-tau*eye(n))'\b */
708:     for (i=0;i<n;i++)
709:       B[i+i*ld] -= tau;
710:     LAPACKgetrf_(&n,&n,B,&ld,ipiv,&info);
711:     if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
712:     if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
713:     PetscLogFlops(2.0*n*n*n/3.0);
714:     LAPACKgetrs_("C",&n,&one,B,&ld,ipiv,g,&ld,&info);
715:     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
716:     PetscLogFlops(2.0*n*n-n);

718:     /* A = A + g*b' */
719:     for (i=0;i<n;i++)
720:       A[i+(n-1)*ld] += g[i]*beta;

722:   } else { /* recover */

725:     DSAllocateWork_Private(ds,ld,0,0);
726:     ghat = ds->work;
727:     Q    = ds->mat[DS_MAT_Q];

729:     /* g^ = -Q(:,idx)'*g */
730:     ncol = PetscBLASIntCast(ds->l+ds->k);
731:     BLASgemv_("C",&n,&ncol,&dmone,Q,&ld,g,&one,&dzero,ghat,&one);

733:     /* A = A + g^*b' */
734:     for (i=0;i<ds->l+ds->k;i++)
735:       for (j=ds->l;j<ds->l+ds->k;j++)
736:         A[i+j*ld] += ghat[i]*Q[n-1+j*ld]*beta;

738:     /* g~ = (I-Q(:,idx)*Q(:,idx)')*g = g+Q(:,idx)*g^ */
739:     BLASgemv_("N",&n,&ncol,&done,Q,&ld,ghat,&one,&done,g,&one);
740:   }

742:   /* Compute gamma factor */
743:   if (gamma) {
744:     gnorm = 0.0;
745:     for (i=0;i<n;i++)
746:       gnorm = gnorm + PetscRealPart(g[i]*PetscConj(g[i]));
747:     *gamma = PetscSqrtReal(1.0+gnorm);
748:   }
749:   return(0);
750: #endif
751: }

753: EXTERN_C_BEGIN
756: PetscErrorCode DSCreate_NHEP(DS ds)
757: {
759:   ds->ops->allocate      = DSAllocate_NHEP;
760:   ds->ops->view          = DSView_NHEP;
761:   ds->ops->vectors       = DSVectors_NHEP;
762:   ds->ops->solve[0]      = DSSolve_NHEP;
763:   ds->ops->sort          = DSSort_NHEP;
764:   ds->ops->truncate      = DSTruncate_NHEP;
765:   ds->ops->update        = DSUpdateExtraRow_NHEP;
766:   ds->ops->cond          = DSCond_NHEP;
767:   ds->ops->transharm     = DSTranslateHarmonic_NHEP;
768:   ds->ops->normalize     = DSNormalize_NHEP;
769:   return(0);
770: }
771: EXTERN_C_END