Actual source code: arnoldi.c

  1: /*                       

  3:    SLEPc eigensolver: "arnoldi"

  5:    Method: Explicitly Restarted Arnoldi

  7:    Algorithm:

  9:        Arnoldi method with explicit restart and deflation.

 11:    References:

 13:        [1] "Arnoldi Methods in SLEPc", SLEPc Technical Report STR-4, 
 14:            available at http://www.grycap.upv.es/slepc.

 16:    Last update: Feb 2009

 18:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 19:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 20:    Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain

 22:    This file is part of SLEPc.
 23:       
 24:    SLEPc is free software: you can redistribute it and/or modify it under  the
 25:    terms of version 3 of the GNU Lesser General Public License as published by
 26:    the Free Software Foundation.

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

 33:    You  should have received a copy of the GNU Lesser General  Public  License
 34:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 35:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 36: */

 38:  #include private/epsimpl.h
 39:  #include slepcblaslapack.h

 41: typedef struct {
 42:   PetscTruth delayed;
 43: } EPS_ARNOLDI;

 47: PetscErrorCode EPSSetUp_ARNOLDI(EPS eps)
 48: {
 50:   PetscInt       N;

 53:   VecGetSize(eps->vec_initial,&N);
 54:   if (eps->ncv) { /* ncv set */
 55:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 56:   }
 57:   else if (eps->mpd) { /* mpd set */
 58:     eps->ncv = PetscMin(N,eps->nev+eps->mpd);
 59:   }
 60:   else { /* neither set: defaults depend on nev being small or large */
 61:     if (eps->nev<500) eps->ncv = PetscMin(N,PetscMax(2*eps->nev,eps->nev+15));
 62:     else { eps->mpd = 500; eps->ncv = PetscMin(N,eps->nev+eps->mpd); }
 63:   }
 64:   if (!eps->mpd) eps->mpd = eps->ncv;
 65:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(1,"The value of ncv must not be larger than nev+mpd");
 66:   if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);
 67:   if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
 68:     SETERRQ(1,"Wrong value of eps->which");

 70:   if (!eps->extraction) {
 71:     EPSSetExtraction(eps,EPS_RITZ);
 72:   }

 74:   EPSAllocateSolution(eps);
 75:   PetscFree(eps->T);
 76:   PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
 77:   if (eps->solverclass==EPS_TWO_SIDE) {
 78:     PetscFree(eps->Tl);
 79:     PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->Tl);
 80:     PetscInfo(eps,"Warning: parameter mpd ignored\n");
 81:   }
 82:   EPSDefaultGetWork(eps,2);
 83:   return(0);
 84: }

 88: /*
 89:    EPSBasicArnoldi - Computes an m-step Arnoldi factorization. The first k
 90:    columns are assumed to be locked and therefore they are not modified. On
 91:    exit, the following relation is satisfied:

 93:                     OP * V - V * H = f * e_m^T

 95:    where the columns of V are the Arnoldi vectors (which are B-orthonormal),
 96:    H is an upper Hessenberg matrix, f is the residual vector and e_m is
 97:    the m-th vector of the canonical basis. The vector f is B-orthogonal to
 98:    the columns of V. On exit, beta contains the B-norm of f and the next 
 99:    Arnoldi vector can be computed as v_{m+1} = f / beta. 
100: */
101: PetscErrorCode EPSBasicArnoldi(EPS eps,PetscTruth trans,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscTruth *breakdown)
102: {
104:   PetscInt       i,j,m = *M;
105:   PetscReal      norm;
106:   PetscScalar    *swork = PETSC_NULL,*hwork = PETSC_NULL;

109:   if (eps->nds+m > 100) { PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&swork); }
110:   if (eps->nds > 0) { PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&hwork); }
111: 
112:   for (j=k;j<m-1;j++) {
113:     if (trans) { STApplyTranspose(eps->OP,V[j],V[j+1]); }
114:     else { STApply(eps->OP,V[j],V[j+1]); }
115:     if (eps->nds > 0) {
116:       IPOrthogonalize(eps->ip,eps->nds+j+1,PETSC_NULL,eps->DSV,V[j+1],hwork,&norm,breakdown,eps->work[0],swork);
117:       for (i=0;i<=j;i++)
118:         H[ldh*j+i] = hwork[eps->nds+i];
119:     } else {
120:       IPOrthogonalize(eps->ip,j+1,PETSC_NULL,V,V[j+1],H+ldh*j,&norm,breakdown,eps->work[0],swork);
121:     }
122:     H[j+1+ldh*j] = norm;
123:     if (*breakdown) {
124:       *M = j+1;
125:       *beta = norm;
126:       if (swork) { PetscFree(swork); }
127:       if (hwork) { PetscFree(hwork); }
128:       return(0);
129:     } else {
130:       VecScale(V[j+1],1/norm);
131:     }
132:   }
133:   if (trans) { STApplyTranspose(eps->OP,V[m-1],f); }
134:   else { STApply(eps->OP,V[m-1],f); }
135:   IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL,eps->work[0],swork);
136:   IPOrthogonalize(eps->ip,m,PETSC_NULL,V,f,H+ldh*(m-1),beta,PETSC_NULL,eps->work[0],swork);
137: 
138:   if (swork) { PetscFree(swork); }
139:   if (hwork) { PetscFree(hwork); }
140:   return(0);
141: }

145: /*
146:    EPSDelayedArnoldi - This function is equivalent to EPSBasicArnoldi but
147:    performs the computation in a different way. The main idea is that
148:    reorthogonalization is delayed to the next Arnoldi step. This version is
149:    more scalable but in some cases convergence may stagnate.
150: */
151: PetscErrorCode EPSDelayedArnoldi(EPS eps,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscTruth *breakdown)
152: {
154:   PetscInt       i,j,m=*M;
155:   Vec            w,u,t;
156:   PetscScalar    shh[100],*lhh,dot,dot2;
157:   PetscReal      norm1=0.0,norm2;

160:   if (m<=100) lhh = shh;
161:   else { PetscMalloc(m*sizeof(PetscScalar),&lhh); }
162:   VecDuplicate(f,&w);
163:   VecDuplicate(f,&u);
164:   VecDuplicate(f,&t);

166:   for (j=k;j<m;j++) {
167:     STApply(eps->OP,V[j],f);
168:     IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL,eps->work[0],PETSC_NULL);

170:     IPMInnerProductBegin(eps->ip,f,j+1,V,H+ldh*j);
171:     if (j>k) {
172:       IPMInnerProductBegin(eps->ip,V[j],j,V,lhh);
173:       IPInnerProductBegin(eps->ip,V[j],V[j],&dot);
174:     }
175:     if (j>k+1) {
176:       IPNormBegin(eps->ip,u,&norm2);
177:       VecDotBegin(u,V[j-2],&dot2);
178:     }
179: 
180:     IPMInnerProductEnd(eps->ip,f,j+1,V,H+ldh*j);
181:     if (j>k) {
182:       IPMInnerProductEnd(eps->ip,V[j],j,V,lhh);
183:       IPInnerProductEnd(eps->ip,V[j],V[j],&dot);
184:     }
185:     if (j>k+1) {
186:       IPNormEnd(eps->ip,u,&norm2);
187:       VecDotEnd(u,V[j-2],&dot2);
188:       if (PetscAbsScalar(dot2/norm2) > PETSC_MACHINE_EPSILON) {
189:         *breakdown = PETSC_TRUE;
190:         *M = j-1;
191:         *beta = norm2;

193:         if (m>100) { PetscFree(lhh); }
194:         VecDestroy(w);
195:         VecDestroy(u);
196:         VecDestroy(t);
197:         return(0);
198:       }
199:     }
200: 
201:     if (j>k) {
202:       norm1 = sqrt(PetscRealPart(dot));
203:       for (i=0;i<j;i++)
204:         H[ldh*j+i] = H[ldh*j+i]/norm1;
205:       H[ldh*j+j] = H[ldh*j+j]/dot;
206: 
207:       VecCopy(V[j],t);
208:       VecScale(V[j],1.0/norm1);
209:       VecScale(f,1.0/norm1);
210:     }

212:     VecSet(w,0.0);
213:     VecMAXPY(w,j+1,H+ldh*j,V);
214:     VecAXPY(f,-1.0,w);

216:     if (j>k) {
217:       VecSet(w,0.0);
218:       VecMAXPY(w,j,lhh,V);
219:       VecAXPY(t,-1.0,w);
220:       for (i=0;i<j;i++)
221:         H[ldh*(j-1)+i] += lhh[i];
222:     }

224:     if (j>k+1) {
225:       VecCopy(u,V[j-1]);
226:       VecScale(V[j-1],1.0/norm2);
227:       H[ldh*(j-2)+j-1] = norm2;
228:     }

230:     if (j<m-1) {
231:       VecCopy(f,V[j+1]);
232:       VecCopy(t,u);
233:     }
234:   }

236:   IPNorm(eps->ip,t,&norm2);
237:   VecScale(t,1.0/norm2);
238:   VecCopy(t,V[m-1]);
239:   H[ldh*(m-2)+m-1] = norm2;

241:   IPMInnerProduct(eps->ip,f,m,V,lhh);
242: 
243:   VecSet(w,0.0);
244:   VecMAXPY(w,m,lhh,V);
245:   VecAXPY(f,-1.0,w);
246:   for (i=0;i<m;i++)
247:     H[ldh*(m-1)+i] += lhh[i];

249:   IPNorm(eps->ip,f,beta);
250:   VecScale(f,1.0 / *beta);
251:   *breakdown = PETSC_FALSE;
252: 
253:   if (m>100) { PetscFree(lhh); }
254:   VecDestroy(w);
255:   VecDestroy(u);
256:   VecDestroy(t);

258:   return(0);
259: }

263: /*
264:    EPSDelayedArnoldi1 - This function is similar to EPSDelayedArnoldi1,
265:    but without reorthogonalization (only delayed normalization).
266: */
267: PetscErrorCode EPSDelayedArnoldi1(EPS eps,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscTruth *breakdown)
268: {
270:   PetscInt       i,j,m=*M;
271:   Vec            w;
272:   PetscScalar    dot;
273:   PetscReal      norm=0.0;

276:   VecDuplicate(f,&w);

278:   for (j=k;j<m;j++) {
279:     STApply(eps->OP,V[j],f);
280:     IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL,eps->work[0],PETSC_NULL);

282:     IPMInnerProductBegin(eps->ip,f,j+1,V,H+ldh*j);
283:     if (j>k) {
284:       IPInnerProductBegin(eps->ip,V[j],V[j],&dot);
285:     }
286: 
287:     IPMInnerProductEnd(eps->ip,f,j+1,V,H+ldh*j);
288:     if (j>k) {
289:       IPInnerProductEnd(eps->ip,V[j],V[j],&dot);
290:     }
291: 
292:     if (j>k) {
293:       norm = sqrt(PetscRealPart(dot));
294:       VecScale(V[j],1.0/norm);
295:       H[ldh*(j-1)+j] = norm;

297:       for (i=0;i<j;i++)
298:         H[ldh*j+i] = H[ldh*j+i]/norm;
299:       H[ldh*j+j] = H[ldh*j+j]/dot;
300:       VecScale(f,1.0/norm);
301:     }

303:     VecSet(w,0.0);
304:     VecMAXPY(w,j+1,H+ldh*j,V);
305:     VecAXPY(f,-1.0,w);

307:     if (j<m-1) {
308:       VecCopy(f,V[j+1]);
309:     }
310:   }

312:   IPNorm(eps->ip,f,beta);
313:   VecScale(f,1.0 / *beta);
314:   *breakdown = PETSC_FALSE;
315: 
316:   VecDestroy(w);
317:   return(0);
318: }

322: /*
323:    EPSArnoldiResiduals - Computes the 2-norm of the residual vectors from
324:    the information provided by the m-step Arnoldi factorization,

326:                     OP * V - V * H = f * e_m^T

328:    For the approximate eigenpair (k_i,V*y_i), the residual norm is computed as
329:    |beta*y(end,i)| where beta is the norm of f and y is the corresponding 
330:    eigenvector of H.
331: */
332: PetscErrorCode ArnoldiResiduals(PetscScalar *H,PetscInt ldh_,PetscScalar *U,PetscReal beta,PetscInt nconv,PetscInt ncv_,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscScalar *work)
333: {
334: #if defined(SLEPC_MISSING_LAPACK_TREVC)
336:   SETERRQ(PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable.");
337: #else
339:   PetscInt       i;
340:   PetscBLASInt   mout,info,ldh,ncv;
341:   PetscScalar    *Y=work+4*ncv_;
342:   PetscReal      w;
343: #if defined(PETSC_USE_COMPLEX)
344:   PetscReal      *rwork=(PetscReal*)(work+3*ncv_);
345: #endif

348:   ldh = PetscBLASIntCast(ldh_);
349:   ncv = PetscBLASIntCast(ncv_);

351:   /* Compute eigenvectors Y of H */
352:   PetscMemcpy(Y,U,ncv*ncv*sizeof(PetscScalar));
353:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
354: #if !defined(PETSC_USE_COMPLEX)
355:   LAPACKtrevc_("R","B",PETSC_NULL,&ncv,H,&ldh,PETSC_NULL,&ncv,Y,&ncv,&ncv,&mout,work,&info);
356: #else
357:   LAPACKtrevc_("R","B",PETSC_NULL,&ncv,H,&ldh,PETSC_NULL,&ncv,Y,&ncv,&ncv,&mout,work,rwork,&info);
358: #endif
359:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
360:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);

362:   /* Compute residual norm estimates as beta*abs(Y(m,:)) */
363:   for (i=nconv;i<ncv;i++) {
364: #if !defined(PETSC_USE_COMPLEX)
365:     if (eigi[i] != 0 && i<ncv-1) {
366:       errest[i] = beta*SlepcAbsEigenvalue(Y[i*ncv+ncv-1],Y[(i+1)*ncv+ncv-1]);
367:       w = SlepcAbsEigenvalue(eigr[i],eigi[i]);
368:       if (w > errest[i])
369:         errest[i] = errest[i] / w;
370:       errest[i+1] = errest[i];
371:       i++;
372:     } else
373: #endif
374:     {
375:       errest[i] = beta*PetscAbsScalar(Y[i*ncv+ncv-1]);
376:       w = PetscAbsScalar(eigr[i]);
377:       if (w > errest[i])
378:         errest[i] = errest[i] / w;
379:     }
380:   }
381:   return(0);
382: #endif
383: }

387: /*
388:    EPSProjectedArnoldi - Solves the projected eigenproblem.

390:    On input:
391:      S is the projected matrix (leading dimension is lds)

393:    On output:
394:      S has (real) Schur form with diagonal blocks sorted appropriately
395:      Q contains the corresponding Schur vectors (order n, leading dimension n)
396: */
397: PetscErrorCode EPSProjectedArnoldi(EPS eps,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt n)
398: {
400:   PetscInt       i;

403:   /* Initialize orthogonal matrix */
404:   PetscMemzero(Q,n*n*sizeof(PetscScalar));
405:   for (i=0;i<n;i++)
406:     Q[i*(n+1)] = 1.0;
407:   /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
408:   EPSDenseSchur(n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);
409:   /* Sort the remaining columns of the Schur form */
410:   if (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC) {
411:     EPSSortDenseSchurTarget(n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi,eps->target,eps->which);
412:   } else {
413:     EPSSortDenseSchur(n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi,eps->which);
414:   }
415:   return(0);
416: }

420: /*
421:    EPSUpdateVectors - Computes approximate Schur vectors (or eigenvectors) by
422:    either Ritz extraction (U=U*Q) or refined Ritz extraction 

424:    On input:
425:      n is the size of U
426:      U is the orthogonal basis of the subspace used for projecting
427:      s is the index of the first vector computed
428:      e+1 is the index of the last vector computed
429:      Q contains the corresponding Schur vectors of the projected matrix (size n x n, leading dimension ldq)
430:      H is the (extended) projected matrix (size n+1 x n, leading dimension ldh)

432:    On output:
433:      v is the resulting vector
434: */
435: PetscErrorCode EPSUpdateVectors(EPS eps,PetscInt n_,Vec *U,PetscInt s,PetscInt e,PetscScalar *Q,PetscInt ldq,PetscScalar *H,PetscInt ldh_)
436: {
437: #if defined(PETSC_MISSING_LAPACK_GESVD) 
438:   SETERRQ(PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable.");
439: #else
441:   PetscTruth     isrefined;
442:   PetscInt       i,j,k;
443:   PetscBLASInt   n1,lwork,idummy=1,info,n=n_,ldh=ldh_;
444:   PetscScalar    *B,sdummy,*work;
445:   PetscReal      *sigma;

448:   isrefined = (eps->extraction==EPS_REFINED || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
449:   if (isrefined) {
450:     /* Refined Ritz extraction */
451:     n1 = n+1;
452:     PetscMalloc(n1*n*sizeof(PetscScalar),&B);
453:     PetscMalloc(6*n*sizeof(PetscReal),&sigma);
454:     lwork = 10*n;
455:     PetscMalloc(lwork*sizeof(PetscScalar),&work);
456: 
457:     for (k=s;k<e;k++) {
458:       /* copy H to B */
459:       for (i=0;i<=n;i++) {
460:         for (j=0;j<n;j++) {
461:           B[i+j*n1] = H[i+j*ldh];
462:         }
463:       }
464:       /* subtract ritz value from diagonal of B^ */
465:       for (i=0;i<n;i++) {
466:         B[i+i*n1] -= eps->eigr[k];  /* MISSING: complex case */
467:       }
468:       /* compute SVD of [H-mu*I] */
469:   #if !defined(PETSC_USE_COMPLEX)
470:       LAPACKgesvd_("N","O",&n1,&n,B,&n1,sigma,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,&info);
471:   #else
472:       LAPACKgesvd_("N","O",&n1,&n,B,&n1,sigma,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,sigma+n,&info);
473:   #endif
474:       if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGESVD %d",info);
475:       /* the smallest singular value is the new error estimate */
476:       eps->errest[k] = sigma[n-1];
477:       /* update vector with right singular vector associated to smallest singular value */
478:       for (i=0;i<n;i++)
479:         Q[k*ldq+i] = B[n-1+i*n1];
480:     }
481:     /* free workspace */
482:     PetscFree(B);
483:     PetscFree(sigma);
484:     PetscFree(work);
485:   }
486:   /* Ritz extraction: v = U*q */
487:   SlepcUpdateVectors(n_,U,s,e,Q,ldq,PETSC_FALSE);
488:   return(0);
489: #endif
490: }

494: PetscErrorCode EPSSolve_ARNOLDI(EPS eps)
495: {
497:   PetscInt       i,k,nv;
498:   Vec            f=eps->work[1];
499:   PetscScalar    *H=eps->T,*U,*g,*work,*Hcopy;
500:   PetscReal      beta,gnorm;
501:   PetscTruth     breakdown;
502:   IPOrthogonalizationRefinementType orthog_ref;
503:   EPS_ARNOLDI    *arnoldi = (EPS_ARNOLDI *)eps->data;

506:   PetscMemzero(eps->T,eps->ncv*eps->ncv*sizeof(PetscScalar));
507:   PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&U);
508:   PetscMalloc((eps->ncv+4)*eps->ncv*sizeof(PetscScalar),&work);
509:   if (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC) {
510:     PetscMalloc(eps->ncv*sizeof(PetscScalar),&g);
511:   }
512:   if (eps->extraction==EPS_REFINED || eps->extraction==EPS_REFINED_HARMONIC) {
513:     PetscMalloc((eps->ncv+1)*eps->ncv*sizeof(PetscScalar),&Hcopy);
514:   }
515: 
516:   IPGetOrthogonalization(eps->ip,PETSC_NULL,&orthog_ref,PETSC_NULL);

518:   /* Get the starting Arnoldi vector */
519:   EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
520: 
521:   /* Restart loop */
522:   while (eps->reason == EPS_CONVERGED_ITERATING) {
523:     eps->its++;

525:     /* Compute an nv-step Arnoldi factorization */
526:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
527:     if (!arnoldi->delayed) {
528:       EPSBasicArnoldi(eps,PETSC_FALSE,H,eps->ncv,eps->V,eps->nconv,&nv,f,&beta,&breakdown);
529:     } else if (orthog_ref == IP_ORTH_REFINE_NEVER) {
530:       EPSDelayedArnoldi1(eps,H,eps->ncv,eps->V,eps->nconv,&nv,f,&beta,&breakdown);
531:     } else {
532:       EPSDelayedArnoldi(eps,H,eps->ncv,eps->V,eps->nconv,&nv,f,&beta,&breakdown);
533:     }

535:     if (eps->extraction==EPS_REFINED || eps->extraction==EPS_REFINED_HARMONIC) {
536:       PetscMemcpy(Hcopy,H,eps->ncv*eps->ncv*sizeof(PetscScalar));
537:       for (i=0;i<nv-1;i++) Hcopy[nv+i*eps->ncv] = 0.0;
538:       Hcopy[nv+(nv-1)*eps->ncv] = beta;
539:     }

541:     /* Compute translation of Krylov decomposition if harmonic extraction used */
542:     if (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC) {
543:       EPSTranslateHarmonic(nv,H,eps->ncv,eps->target,(PetscScalar)beta,g,work);
544:     }

546:     /* Solve projected problem and compute residual norm estimates */
547:     EPSProjectedArnoldi(eps,H,eps->ncv,U,nv);
548:     ArnoldiResiduals(H,eps->ncv,U,beta,eps->nconv,nv,eps->eigr,eps->eigi,eps->errest,work);
549: 
550:     /* Fix residual norms if harmonic */
551:     if (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC) {
552:       gnorm = 0.0;
553:       for (i=0;i<nv;i++)
554:         gnorm = gnorm + PetscRealPart(g[i]*PetscConj(g[i]));
555:       for (i=eps->nconv;i<nv;i++)
556:         eps->errest[i] *= sqrt(1.0+gnorm);
557:     }

559:     /* Lock converged eigenpairs and update the corresponding vectors,
560:        including the restart vector: V(:,idx) = V*U(:,idx) */
561:     k = eps->nconv;
562:     while (k<nv && eps->errest[k]<eps->tol) k++;
563:     EPSUpdateVectors(eps,nv,eps->V,eps->nconv,PetscMin(k+1,nv),U,nv,Hcopy,eps->ncv);
564:     eps->nconv = k;

566:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
567:     if (breakdown) {
568:       PetscInfo2(eps,"Breakdown in Arnoldi method (it=%i norm=%g)\n",eps->its,beta);
569:       EPSGetStartVector(eps,k,eps->V[k],&breakdown);
570:       if (breakdown) {
571:         eps->reason = EPS_DIVERGED_BREAKDOWN;
572:         PetscInfo(eps,"Unable to generate more start vectors\n");
573:       }
574:     }
575:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
576:     if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
577:   }
578: 
579:   PetscFree(U);
580:   PetscFree(work);
581:   if (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC) {
582:     PetscFree(g);
583:   }
584:   if (eps->extraction==EPS_REFINED || eps->extraction==EPS_REFINED_HARMONIC) {
585:     PetscFree(Hcopy);
586:   }
587:   return(0);
588: }

592: PetscErrorCode EPSSetFromOptions_ARNOLDI(EPS eps)
593: {
595:   EPS_ARNOLDI    *arnoldi = (EPS_ARNOLDI *)eps->data;

598:   PetscOptionsHead("ARNOLDI options");
599:   PetscOptionsTruth("-eps_arnoldi_delayed","Arnoldi with delayed reorthogonalization","EPSArnoldiSetDelayed",PETSC_FALSE,&arnoldi->delayed,PETSC_NULL);
600:   PetscOptionsTail();
601:   return(0);
602: }

607: PetscErrorCode EPSArnoldiSetDelayed_ARNOLDI(EPS eps,PetscTruth delayed)
608: {
609:   EPS_ARNOLDI    *arnoldi = (EPS_ARNOLDI *)eps->data;

612:   arnoldi->delayed = delayed;
613:   return(0);
614: }

619: /*@
620:    EPSArnoldiSetDelayed - Activates or deactivates delayed reorthogonalization 
621:    in the Arnoldi iteration. 

623:    Collective on EPS

625:    Input Parameters:
626: +  eps - the eigenproblem solver context
627: -  delayed - boolean flag

629:    Options Database Key:
630: .  -eps_arnoldi_delayed - Activates delayed reorthogonalization in Arnoldi
631:    
632:    Note:
633:    Delayed reorthogonalization is an aggressive optimization for the Arnoldi
634:    eigensolver than may provide better scalability, but sometimes makes the
635:    solver converge less than the default algorithm.

637:    Level: advanced

639: .seealso: EPSArnoldiGetDelayed()
640: @*/
641: PetscErrorCode EPSArnoldiSetDelayed(EPS eps,PetscTruth delayed)
642: {
643:   PetscErrorCode ierr, (*f)(EPS,PetscTruth);

647:   PetscObjectQueryFunction((PetscObject)eps,"EPSArnoldiSetDelayed_C",(void (**)())&f);
648:   if (f) {
649:     (*f)(eps,delayed);
650:   }
651:   return(0);
652: }

657: PetscErrorCode EPSArnoldiGetDelayed_ARNOLDI(EPS eps,PetscTruth *delayed)
658: {
659:   EPS_ARNOLDI    *arnoldi = (EPS_ARNOLDI *)eps->data;

662:   *delayed = arnoldi->delayed;
663:   return(0);
664: }

669: /*@C
670:    EPSArnoldiGetDelayed - Gets the type of reorthogonalization used during the Arnoldi
671:    iteration. 

673:    Collective on EPS

675:    Input Parameter:
676: .  eps - the eigenproblem solver context

678:    Input Parameter:
679: .  delayed - boolean flag indicating if delayed reorthogonalization has been enabled

681:    Level: advanced

683: .seealso: EPSArnoldiSetDelayed()
684: @*/
685: PetscErrorCode EPSArnoldiGetDelayed(EPS eps,PetscTruth *delayed)
686: {
687:   PetscErrorCode ierr, (*f)(EPS,PetscTruth*);

691:   PetscObjectQueryFunction((PetscObject)eps,"EPSArnoldiGetDelayed_C",(void (**)())&f);
692:   if (f) {
693:     (*f)(eps,delayed);
694:   }
695:   return(0);
696: }

700: PetscErrorCode EPSView_ARNOLDI(EPS eps,PetscViewer viewer)
701: {
703:   PetscTruth     isascii;
704:   EPS_ARNOLDI    *arnoldi = (EPS_ARNOLDI *)eps->data;

707:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
708:   if (!isascii) {
709:     SETERRQ1(1,"Viewer type %s not supported for EPSARNOLDI",((PetscObject)viewer)->type_name);
710:   }
711:   if (arnoldi->delayed) {
712:     PetscViewerASCIIPrintf(viewer,"using delayed reorthogonalization\n");
713:   }
714:   return(0);
715: }

717: EXTERN PetscErrorCode EPSSolve_TS_ARNOLDI(EPS);

722: PetscErrorCode EPSCreate_ARNOLDI(EPS eps)
723: {
725:   EPS_ARNOLDI    *arnoldi;
726: 
728:   PetscNew(EPS_ARNOLDI,&arnoldi);
729:   PetscLogObjectMemory(eps,sizeof(EPS_ARNOLDI));
730:   eps->data                      = (void *)arnoldi;
731:   eps->ops->solve                = EPSSolve_ARNOLDI;
732:   eps->ops->solvets              = EPSSolve_TS_ARNOLDI;
733:   eps->ops->setup                = EPSSetUp_ARNOLDI;
734:   eps->ops->setfromoptions       = EPSSetFromOptions_ARNOLDI;
735:   eps->ops->destroy              = EPSDestroy_Default;
736:   eps->ops->view                 = EPSView_ARNOLDI;
737:   eps->ops->backtransform        = EPSBackTransform_Default;
738:   eps->ops->computevectors       = EPSComputeVectors_Schur;
739:   arnoldi->delayed               = PETSC_FALSE;
740:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSArnoldiSetDelayed_C","EPSArnoldiSetDelayed_ARNOLDI",EPSArnoldiSetDelayed_ARNOLDI);
741:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSArnoldiGetDelayed_C","EPSArnoldiGetDelayed_ARNOLDI",EPSArnoldiGetDelayed_ARNOLDI);
742:   return(0);
743: }