Actual source code: arnoldi.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  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:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 17:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 18:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

 20:    This file is part of SLEPc.

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

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

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

 36: #include <slepc-private/epsimpl.h>                /*I "slepceps.h" I*/

 38: PetscErrorCode EPSSolve_Arnoldi(EPS);

 40: typedef struct {
 41:   PetscBool delayed;
 42: } EPS_ARNOLDI;

 46: PetscErrorCode EPSSetUp_Arnoldi(EPS eps)
 47: {

 51:   EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
 52:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
 53:   if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 54:   if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
 55:   if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");

 57:   if (!eps->extraction) {
 58:     EPSSetExtraction(eps,EPS_RITZ);
 59:   }
 60:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");

 62:   EPSAllocateSolution(eps,1);
 63:   EPS_SetInnerProduct(eps);
 64:   DSSetType(eps->ds,DSNHEP);
 65:   if (eps->extraction==EPS_REFINED || eps->extraction==EPS_REFINED_HARMONIC) {
 66:     DSSetRefined(eps->ds,PETSC_TRUE);
 67:   }
 68:   DSSetExtraRow(eps->ds,PETSC_TRUE);
 69:   DSAllocate(eps->ds,eps->ncv+1);

 71:   /* dispatch solve method */
 72:   if (eps->isgeneralized && eps->ishermitian && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method does not work for indefinite problems");
 73:   eps->ops->solve = EPSSolve_Arnoldi;
 74:   return(0);
 75: }

 77: #if 0
 80: /*
 81:    EPSDelayedArnoldi - This function is equivalent to EPSBasicArnoldi but
 82:    performs the computation in a different way. The main idea is that
 83:    reorthogonalization is delayed to the next Arnoldi step. This version is
 84:    more scalable but in some cases convergence may stagnate.
 85: */
 86: PetscErrorCode EPSDelayedArnoldi(EPS eps,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscBool *breakdown)
 87: {
 89:   PetscInt       i,j,m=*M;
 90:   Vec            u,t;
 91:   PetscScalar    shh[100],*lhh,dot,dot2;
 92:   PetscReal      norm1=0.0,norm2;

 95:   if (m<=100) lhh = shh;
 96:   else {
 97:     PetscMalloc1(m,&lhh);
 98:   }
 99:   VecDuplicate(f,&u);
100:   VecDuplicate(f,&t);

102:   for (j=k;j<m;j++) {
103:     STApply(eps->st,V[j],f);
104:     IPOrthogonalize(eps->ip,0,NULL,eps->nds,NULL,eps->defl,f,NULL,NULL,NULL);

106:     IPMInnerProductBegin(eps->ip,f,j+1,V,H+ldh*j);
107:     if (j>k) {
108:       IPMInnerProductBegin(eps->ip,V[j],j,V,lhh);
109:       IPInnerProductBegin(eps->ip,V[j],V[j],&dot);
110:     }
111:     if (j>k+1) {
112:       IPNormBegin(eps->ip,u,&norm2);
113:       VecDotBegin(u,V[j-2],&dot2);
114:     }

116:     IPMInnerProductEnd(eps->ip,f,j+1,V,H+ldh*j);
117:     if (j>k) {
118:       IPMInnerProductEnd(eps->ip,V[j],j,V,lhh);
119:       IPInnerProductEnd(eps->ip,V[j],V[j],&dot);
120:     }
121:     if (j>k+1) {
122:       IPNormEnd(eps->ip,u,&norm2);
123:       VecDotEnd(u,V[j-2],&dot2);
124:       if (PetscAbsScalar(dot2/norm2) > PETSC_MACHINE_EPSILON) {
125:         *breakdown = PETSC_TRUE;
126:         *M = j-1;
127:         *beta = norm2;

129:         if (m>100) { PetscFree(lhh); }
130:         VecDestroy(&u);
131:         VecDestroy(&t);
132:         return(0);
133:       }
134:     }

136:     if (j>k) {
137:       norm1 = PetscSqrtReal(PetscRealPart(dot));
138:       for (i=0;i<j;i++)
139:         H[ldh*j+i] = H[ldh*j+i]/norm1;
140:       H[ldh*j+j] = H[ldh*j+j]/dot;

142:       VecCopy(V[j],t);
143:       VecScale(V[j],1.0/norm1);
144:       VecScale(f,1.0/norm1);
145:     }

147:     SlepcVecMAXPBY(f,1.0,-1.0,j+1,H+ldh*j,V);

149:     if (j>k) {
150:       SlepcVecMAXPBY(t,1.0,-1.0,j,lhh,V);
151:       for (i=0;i<j;i++)
152:         H[ldh*(j-1)+i] += lhh[i];
153:     }

155:     if (j>k+1) {
156:       VecCopy(u,V[j-1]);
157:       VecScale(V[j-1],1.0/norm2);
158:       H[ldh*(j-2)+j-1] = norm2;
159:     }

161:     if (j<m-1) {
162:       VecCopy(f,V[j+1]);
163:       VecCopy(t,u);
164:     }
165:   }

167:   IPNorm(eps->ip,t,&norm2);
168:   VecScale(t,1.0/norm2);
169:   VecCopy(t,V[m-1]);
170:   H[ldh*(m-2)+m-1] = norm2;

172:   IPMInnerProduct(eps->ip,f,m,V,lhh);

174:   SlepcVecMAXPBY(f,1.0,-1.0,m,lhh,V);
175:   for (i=0;i<m;i++)
176:     H[ldh*(m-1)+i] += lhh[i];

178:   IPNorm(eps->ip,f,beta);
179:   VecScale(f,1.0 / *beta);
180:   *breakdown = PETSC_FALSE;

182:   if (m>100) { PetscFree(lhh); }
183:   VecDestroy(&u);
184:   VecDestroy(&t);
185:   return(0);
186: }

190: /*
191:    EPSDelayedArnoldi1 - This function is similar to EPSDelayedArnoldi,
192:    but without reorthogonalization (only delayed normalization).
193: */
194: PetscErrorCode EPSDelayedArnoldi1(EPS eps,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscBool *breakdown)
195: {
197:   PetscInt       i,j,m=*M;
198:   PetscScalar    dot;
199:   PetscReal      norm=0.0;

202:   for (j=k;j<m;j++) {
203:     STApply(eps->st,V[j],f);
204:     IPOrthogonalize(eps->ip,0,NULL,eps->nds,NULL,eps->defl,f,NULL,NULL,NULL);

206:     IPMInnerProductBegin(eps->ip,f,j+1,V,H+ldh*j);
207:     if (j>k) {
208:       IPInnerProductBegin(eps->ip,V[j],V[j],&dot);
209:     }

211:     IPMInnerProductEnd(eps->ip,f,j+1,V,H+ldh*j);
212:     if (j>k) {
213:       IPInnerProductEnd(eps->ip,V[j],V[j],&dot);
214:     }

216:     if (j>k) {
217:       norm = PetscSqrtReal(PetscRealPart(dot));
218:       VecScale(V[j],1.0/norm);
219:       H[ldh*(j-1)+j] = norm;

221:       for (i=0;i<j;i++)
222:         H[ldh*j+i] = H[ldh*j+i]/norm;
223:       H[ldh*j+j] = H[ldh*j+j]/dot;
224:       VecScale(f,1.0/norm);
225:     }

227:     SlepcVecMAXPBY(f,1.0,-1.0,j+1,H+ldh*j,V);

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

234:   IPNorm(eps->ip,f,beta);
235:   VecScale(f,1.0 / *beta);
236:   *breakdown = PETSC_FALSE;
237:   return(0);
238: }
239: #endif

243: PetscErrorCode EPSSolve_Arnoldi(EPS eps)
244: {
245:   PetscErrorCode     ierr;
246:   PetscInt           k,nv,ld;
247:   Mat                U;
248:   PetscScalar        *H,*X;
249:   PetscReal          beta,gamma=1.0;
250:   PetscBool          breakdown,harmonic,refined;
251:   BVOrthogRefineType orthog_ref;
252:   EPS_ARNOLDI        *arnoldi = (EPS_ARNOLDI*)eps->data;

255:   DSGetLeadingDimension(eps->ds,&ld);
256:   DSGetRefined(eps->ds,&refined);
257:   harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
258:   BVGetOrthogonalization(eps->V,NULL,&orthog_ref,NULL);

260:   /* Get the starting Arnoldi vector */
261:   EPSGetStartVector(eps,0,NULL);

263:   /* Restart loop */
264:   while (eps->reason == EPS_CONVERGED_ITERATING) {
265:     eps->its++;

267:     /* Compute an nv-step Arnoldi factorization */
268:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
269:     DSSetDimensions(eps->ds,nv,0,eps->nconv,0);
270:     DSGetArray(eps->ds,DS_MAT_A,&H);
271:     if (!arnoldi->delayed) {
272:       EPSBasicArnoldi(eps,PETSC_FALSE,H,ld,eps->nconv,&nv,&beta,&breakdown);
273:     } else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Not implemented");
274:     /*if (orthog_ref == BV_ORTHOG_REFINE_NEVER) {
275:       EPSDelayedArnoldi1(eps,H,ld,eps->V,eps->nconv,&nv,f,&beta,&breakdown);
276:     } else {
277:       EPSDelayedArnoldi(eps,H,ld,eps->V,eps->nconv,&nv,f,&beta,&breakdown);
278:     }*/
279:     DSRestoreArray(eps->ds,DS_MAT_A,&H);
280:     DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
281:     BVSetActiveColumns(eps->V,eps->nconv,nv);

283:     /* Compute translation of Krylov decomposition if harmonic extraction used */
284:     if (harmonic) {
285:       DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,NULL,&gamma);
286:     }

288:     /* Solve projected problem */
289:     DSSolve(eps->ds,eps->eigr,eps->eigi);
290:     DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
291:     DSUpdateExtraRow(eps->ds);

293:     /* Check convergence */
294:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,gamma,&k);
295:     if (refined) {
296:       DSGetArray(eps->ds,DS_MAT_X,&X);
297:       BVMultColumn(eps->V,1.0,0.0,k,X+k*ld);
298:       DSRestoreArray(eps->ds,DS_MAT_X,&X);
299:       DSGetMat(eps->ds,DS_MAT_Q,&U);
300:       BVMultInPlace(eps->V,U,eps->nconv,nv);
301:       MatDestroy(&U);
302:       BVOrthogonalizeColumn(eps->V,k,NULL,NULL,NULL);
303:     } else {
304:       DSGetMat(eps->ds,DS_MAT_Q,&U);
305:       BVMultInPlace(eps->V,U,eps->nconv,nv);
306:       MatDestroy(&U);
307:     }
308:     eps->nconv = k;

310:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
311:     if (breakdown && k<eps->nev) {
312:       PetscInfo2(eps,"Breakdown in Arnoldi method (it=%D norm=%g)\n",eps->its,(double)beta);
313:       EPSGetStartVector(eps,k,&breakdown);
314:       if (breakdown) {
315:         eps->reason = EPS_DIVERGED_BREAKDOWN;
316:         PetscInfo(eps,"Unable to generate more start vectors\n");
317:       }
318:     }
319:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
320:     if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
321:   }

323:   /* truncate Schur decomposition and change the state to raw so that
324:      PSVectors() computes eigenvectors from scratch */
325:   DSSetDimensions(eps->ds,eps->nconv,0,0,0);
326:   DSSetState(eps->ds,DS_STATE_RAW);
327:   return(0);
328: }

332: PetscErrorCode EPSSetFromOptions_Arnoldi(EPS eps)
333: {
335:   PetscBool      set,val;
336:   EPS_ARNOLDI    *arnoldi = (EPS_ARNOLDI*)eps->data;

339:   PetscOptionsHead("EPS Arnoldi Options");
340:   PetscOptionsBool("-eps_arnoldi_delayed","Arnoldi with delayed reorthogonalization","EPSArnoldiSetDelayed",arnoldi->delayed,&val,&set);
341:   if (set) {
342:     EPSArnoldiSetDelayed(eps,val);
343:   }
344:   PetscOptionsTail();
345:   return(0);
346: }

350: static PetscErrorCode EPSArnoldiSetDelayed_Arnoldi(EPS eps,PetscBool delayed)
351: {
352:   EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;

355:   arnoldi->delayed = delayed;
356:   return(0);
357: }

361: /*@
362:    EPSArnoldiSetDelayed - Activates or deactivates delayed reorthogonalization
363:    in the Arnoldi iteration.

365:    Logically Collective on EPS

367:    Input Parameters:
368: +  eps - the eigenproblem solver context
369: -  delayed - boolean flag

371:    Options Database Key:
372: .  -eps_arnoldi_delayed - Activates delayed reorthogonalization in Arnoldi

374:    Note:
375:    Delayed reorthogonalization is an aggressive optimization for the Arnoldi
376:    eigensolver than may provide better scalability, but sometimes makes the
377:    solver converge less than the default algorithm.

379:    Level: advanced

381: .seealso: EPSArnoldiGetDelayed()
382: @*/
383: PetscErrorCode EPSArnoldiSetDelayed(EPS eps,PetscBool delayed)
384: {

390:   PetscTryMethod(eps,"EPSArnoldiSetDelayed_C",(EPS,PetscBool),(eps,delayed));
391:   return(0);
392: }

396: static PetscErrorCode EPSArnoldiGetDelayed_Arnoldi(EPS eps,PetscBool *delayed)
397: {
398:   EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;

401:   *delayed = arnoldi->delayed;
402:   return(0);
403: }

407: /*@C
408:    EPSArnoldiGetDelayed - Gets the type of reorthogonalization used during the Arnoldi
409:    iteration.

411:    Not Collective

413:    Input Parameter:
414: .  eps - the eigenproblem solver context

416:    Input Parameter:
417: .  delayed - boolean flag indicating if delayed reorthogonalization has been enabled

419:    Level: advanced

421: .seealso: EPSArnoldiSetDelayed()
422: @*/
423: PetscErrorCode EPSArnoldiGetDelayed(EPS eps,PetscBool *delayed)
424: {

430:   PetscTryMethod(eps,"EPSArnoldiGetDelayed_C",(EPS,PetscBool*),(eps,delayed));
431:   return(0);
432: }

436: PetscErrorCode EPSDestroy_Arnoldi(EPS eps)
437: {

441:   PetscFree(eps->data);
442:   PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiSetDelayed_C",NULL);
443:   PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiGetDelayed_C",NULL);
444:   return(0);
445: }

449: PetscErrorCode EPSView_Arnoldi(EPS eps,PetscViewer viewer)
450: {
452:   PetscBool      isascii;
453:   EPS_ARNOLDI    *arnoldi = (EPS_ARNOLDI*)eps->data;

456:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
457:   if (isascii) {
458:     if (arnoldi->delayed) {
459:       PetscViewerASCIIPrintf(viewer,"  Arnoldi: using delayed reorthogonalization\n");
460:     }
461:   }
462:   return(0);
463: }

467: PETSC_EXTERN PetscErrorCode EPSCreate_Arnoldi(EPS eps)
468: {
469:   EPS_ARNOLDI    *ctx;

473:   PetscNewLog(eps,&ctx);
474:   eps->data = (void*)ctx;

476:   eps->ops->setup                = EPSSetUp_Arnoldi;
477:   eps->ops->setfromoptions       = EPSSetFromOptions_Arnoldi;
478:   eps->ops->destroy              = EPSDestroy_Arnoldi;
479:   eps->ops->view                 = EPSView_Arnoldi;
480:   eps->ops->backtransform        = EPSBackTransform_Default;
481:   eps->ops->computevectors       = EPSComputeVectors_Schur;
482:   PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiSetDelayed_C",EPSArnoldiSetDelayed_Arnoldi);
483:   PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiGetDelayed_C",EPSArnoldiGetDelayed_Arnoldi);
484:   return(0);
485: }