Actual source code: power.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc eigensolver: "power"

 13:    Method: Power Iteration

 15:    Algorithm:

 17:        This solver implements the power iteration for finding dominant
 18:        eigenpairs. It also includes the following well-known methods:
 19:        - Inverse Iteration: when used in combination with shift-and-invert
 20:          spectral transformation.
 21:        - Rayleigh Quotient Iteration (RQI): also with shift-and-invert plus
 22:          a variable shift.

 24:        It can also be used for nonlinear inverse iteration on the problem
 25:        A(x)*x=lambda*B(x)*x, where A and B are not constant but depend on x.

 27:    References:

 29:        [1] "Single Vector Iteration Methods in SLEPc", SLEPc Technical Report
 30:            STR-2, available at http://slepc.upv.es.
 31: */

 33: #include <slepc/private/epsimpl.h>                /*I "slepceps.h" I*/
 34: #include <slepcblaslapack.h>
 35: /* petsc headers */
 36: #include <petscdm.h>
 37: #include <petscsnes.h>

 39: static PetscErrorCode EPSPowerFormFunction_Update(SNES,Vec,Vec,void*);
 40: static PetscErrorCode SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx);

 42: typedef struct {
 43:   EPSPowerShiftType shift_type;
 44:   PetscBool         nonlinear;
 45:   PetscBool         update;
 46:   SNES              snes;
 47:   PetscErrorCode    (*formFunctionB)(SNES,Vec,Vec,void*);
 48:   void              *formFunctionBctx;
 49:   PetscErrorCode    (*formFunctionA)(SNES,Vec,Vec,void*);
 50:   void              *formFunctionActx;
 51: } EPS_POWER;

 53: PetscErrorCode EPSSetUp_Power(EPS eps)
 54: {
 56:   EPS_POWER      *power = (EPS_POWER*)eps->data;
 57:   PetscBool      flg,istrivial;
 58:   STMatMode      mode;
 59:   Mat            A,B;
 60:   Vec            res;
 61:   PetscContainer container;
 62:   PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
 63:   PetscErrorCode (*formJacobianA)(SNES,Vec,Mat,Mat,void*);
 64:   void           *ctx;
 65:   SNESLineSearch linesearch;

 68:   if (eps->ncv) {
 69:     if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
 70:   } else eps->ncv = eps->nev;
 71:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 72:   if (!eps->max_it) eps->max_it = PetscMax(1000*eps->nev,100*eps->n);
 73:   if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
 74:   if (eps->which!=EPS_LARGEST_MAGNITUDE && eps->which !=EPS_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
 75:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
 76:     if (power->nonlinear) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts not allowed in nonlinear problems");
 77:     PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STSINVERT,STCAYLEY,"");
 78:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert or Cayley ST");
 79:     STGetMatMode(eps->st,&mode);
 80:     if (mode == ST_MATMODE_INPLACE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
 81:   }
 82:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
 83:   if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Balancing not supported in this solver");
 84:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
 85:   RGIsTrivial(eps->rg,&istrivial);
 86:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
 87:   EPSAllocateSolution(eps,0);
 88:   EPS_SetInnerProduct(eps);

 90:   if (power->nonlinear) {
 91:     if (eps->nev>1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration cannot compute more than one eigenvalue");
 92:     EPSSetWorkVecs(eps,4);

 94:     /* set up SNES solver */
 95:     if (!power->snes) { EPSPowerGetSNES(eps,&power->snes); }
 96:     EPSGetOperators(eps,&A,&B);

 98:     PetscObjectQueryFunction((PetscObject)A,"formFunction",&formFunctionA);
 99:     if (!formFunctionA) SETERRQ(PetscObjectComm((PetscObject)eps),1,"For nonlinear inverse iteration you must compose a callback function 'formFunction' in matrix A");
100:     PetscObjectQuery((PetscObject)A,"formFunctionCtx",(PetscObject*)&container);
101:     if (container) {
102:       PetscContainerGetPointer(container,&ctx);
103:     } else ctx = NULL;
104:     MatCreateVecs(A,&res,NULL);
105:     power->formFunctionA = formFunctionA;
106:     power->formFunctionActx = ctx;
107:     if (power->update) { SNESSetFunction(power->snes,res,EPSPowerFormFunction_Update,ctx); }
108:     else { SNESSetFunction(power->snes,res,formFunctionA,ctx); }
109:     VecDestroy(&res);

111:     PetscObjectQueryFunction((PetscObject)A,"formJacobian",&formJacobianA);
112:     if (!formJacobianA) SETERRQ(PetscObjectComm((PetscObject)eps),1,"For nonlinear inverse iteration you must compose a callback function 'formJacobian' in matrix A");
113:     PetscObjectQuery((PetscObject)A,"formJacobianCtx",(PetscObject*)&container);
114:     if (container) {
115:       PetscContainerGetPointer(container,&ctx);
116:     } else ctx = NULL;
117:     SNESSetJacobian(power->snes,A,A,formJacobianA,ctx);
118:     SNESSetFromOptions(power->snes);
119:     SNESGetLineSearch(power->snes,&linesearch);
120:     if (power->update) {
121:       SNESLineSearchSetPostCheck(linesearch,SNESLineSearchPostheckFunction,ctx);
122:     }
123:     SNESSetUp(power->snes);
124:     if (B) {
125:       PetscObjectQueryFunction((PetscObject)B,"formFunction",&power->formFunctionB);
126:       PetscObjectQuery((PetscObject)B,"formFunctionCtx",(PetscObject*)&container);
127:       if (power->formFunctionB && container) {
128:         PetscContainerGetPointer(container,&power->formFunctionBctx);
129:       } else power->formFunctionBctx = NULL;
130:     }
131:   } else {
132:     EPSSetWorkVecs(eps,2);
133:     DSSetType(eps->ds,DSNHEP);
134:     DSAllocate(eps->ds,eps->nev);
135:   }
136:   return(0);
137: }


140: /*
141:    Normalize a vector x with respect to a given norm as well as the
142:    sign of the first entry.
143: */
144: static PetscErrorCode Normalize(Vec x,PetscReal norm,PetscScalar *sign)
145: {
146:   PetscErrorCode    ierr;
147:   PetscScalar       alpha = 1.0;
148:   PetscMPIInt       rank;
149:   PetscReal         absv;
150:   const PetscScalar *xx;

153:   MPI_Comm_rank(PetscObjectComm((PetscObject)x),&rank);
154:   if (!rank) {
155:     VecGetArrayRead(x,&xx);
156:     absv = PetscAbsScalar(*xx);
157:     if (absv>10*PETSC_MACHINE_EPSILON) alpha = *xx/absv;
158:     VecRestoreArrayRead(x,&xx);
159:   }
160:   MPI_Bcast(&alpha,1,MPIU_SCALAR,0,PetscObjectComm((PetscObject)x));
161:   if (sign) *sign = alpha;
162:   alpha *= norm;
163:   VecScale(x,1.0/alpha);
164:   return(0);
165: }

167: static PetscErrorCode EPSPowerUpdateFunctionB(EPS eps,Vec x,Vec Bx)
168: {
170:   EPS_POWER      *power = (EPS_POWER*)eps->data;
171:   Mat            B;

174:   STResetMatrixState(eps->st);
175:   EPSGetOperators(eps,NULL,&B);
176:   if (B) {
177:     if (power->formFunctionB) {
178:       (*power->formFunctionB)(power->snes,x,Bx,power->formFunctionBctx);
179:     } else {
180:       MatMult(B,x,Bx);
181:     }
182:   } else {
183:     VecCopy(x,Bx);
184:   }
185:   return(0);
186: }

188: static PetscErrorCode EPSPowerUpdateFunctionA(EPS eps,Vec x,Vec Ax)
189: {
191:   EPS_POWER      *power = (EPS_POWER*)eps->data;
192:   Mat            A;

195:   STResetMatrixState(eps->st);
196:   EPSGetOperators(eps,&A,NULL);
197:   if (A) {
198:     if (power->formFunctionA) {
199:       (*power->formFunctionA)(power->snes,x,Ax,power->formFunctionActx);
200:     } else {
201:       MatMult(A,x,Ax);
202:     }
203:   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NULL,"Matrix A is required for an eigenvalue problem\n");
204:   return(0);
205: }

207: static PetscErrorCode SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx)
208: {
210:   SNES           snes;
211:   EPS            eps;
212:   Vec            oldx;

215:   SNESLineSearchGetSNES(linesearch,&snes);
216:   PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps);
217:   if (!eps) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NULL,"No composed EPS");
218:   oldx = eps->work[3];
219:   VecCopy(x,oldx);
220:   return(0);
221: }

223: static PetscErrorCode EPSPowerFormFunction_Update(SNES snes,Vec x,Vec y,void *ctx)
224: {
226:   EPS            eps;
227:   PetscReal      bx;
228:   Vec            Bx;

231:   PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps);
232:   if (!eps) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NULL,"No composed EPS");
233:   Bx = eps->work[2];
234:   EPSPowerUpdateFunctionB(eps,x,Bx);
235:   VecNorm(Bx,NORM_2,&bx);
236:   Normalize(Bx,bx,NULL);
237:   EPSPowerUpdateFunctionA(eps,x,y);
238:   VecAXPY(y,-1.0,Bx);
239:   return(0);
240: }

242: /*
243:    Use SNES to compute y = A^{-1}*B*x for the nonlinear problem
244: */
245: static PetscErrorCode EPSPowerApply_SNES(EPS eps,Vec x,Vec y)
246: {
248:   EPS_POWER      *power = (EPS_POWER*)eps->data;
249:   Vec            Bx;

252:   VecCopy(x,y);
253:   if (power->update) {
254:     SNESSolve(power->snes,NULL,y);
255:   } else {
256:     Bx = eps->work[2];
257:     SNESSolve(power->snes,Bx,y);
258:   }
259:   return(0);
260: }

262: /*
263:    Use nonlinear inverse power to compute an initial guess
264: */
265: static PetscErrorCode EPSPowerComputeInitialGuess_Update(EPS eps)
266: {
267:   EPS            powereps;
268:   Mat            A,B;
269:   Vec            v1,v2;
270:   SNES           snes;
271:   DM             dm,newdm;

275:   EPSCreate(PetscObjectComm((PetscObject)eps),&powereps);
276:   EPSGetOperators(eps,&A,&B);
277:   EPSSetType(powereps,EPSPOWER);
278:   EPSSetOperators(powereps,A,B);
279:   EPSSetTolerances(powereps,1e-6,4);
280:   EPSSetOptionsPrefix(powereps,((PetscObject)eps)->prefix);
281:   EPSAppendOptionsPrefix(powereps,"init_");
282:   EPSSetProblemType(powereps,EPS_GNHEP);
283:   EPSSetWhichEigenpairs(powereps,EPS_TARGET_MAGNITUDE);
284:   EPSPowerSetNonlinear(powereps,PETSC_TRUE);
285:   STSetType(powereps->st,STSINVERT);
286:   /* attach dm to initial solve */
287:   EPSPowerGetSNES(eps,&snes);
288:   SNESGetDM(snes,&dm);
289:   /* use  dmshell to temporarily store snes context */
290:   DMCreate(PetscObjectComm((PetscObject)eps),&newdm);
291:   DMSetType(newdm,DMSHELL);
292:   DMSetUp(newdm);
293:   DMCopyDMSNES(dm,newdm);
294:   EPSPowerGetSNES(powereps,&snes);
295:   SNESSetDM(snes,dm);
296:   EPSSetFromOptions(powereps);
297:   EPSSolve(powereps);
298:   BVGetColumn(eps->V,0,&v2);
299:   BVGetColumn(powereps->V,0,&v1);
300:   VecCopy(v1,v2);
301:   BVRestoreColumn(powereps->V,0,&v1);
302:   BVRestoreColumn(eps->V,0,&v2);
303:   EPSDestroy(&powereps);
304:   /* restore context back to the old nonlinear solver */
305:   DMCopyDMSNES(newdm,dm);
306:   DMDestroy(&newdm);
307:   return(0);
308: }

310: PetscErrorCode EPSSolve_Power(EPS eps)
311: {
312: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
314:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable");
315: #else
316:   PetscErrorCode      ierr;
317:   EPS_POWER           *power = (EPS_POWER*)eps->data;
318:   PetscInt            k,ld;
319:   Vec                 v,y,e,Bx;
320:   Mat                 A;
321:   KSP                 ksp;
322:   PetscReal           relerr,norm,norm1,rt1,rt2,cs1;
323:   PetscScalar         theta,rho,delta,sigma,alpha2,beta1,sn1,*T,sign;
324:   PetscBool           breakdown;
325:   KSPConvergedReason  reason;
326:   SNESConvergedReason snesreason;

329:   e = eps->work[0];
330:   y = eps->work[1];
331:   if (power->nonlinear) Bx = eps->work[2];
332:   else Bx = NULL;

334:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) { STGetKSP(eps->st,&ksp); }
335:   if (eps->useds) { DSGetLeadingDimension(eps->ds,&ld); }
336:   if (power->nonlinear) {
337:     PetscObjectCompose((PetscObject)power->snes,"eps",(PetscObject)eps);
338:     if (power->update) {
339:       EPSPowerComputeInitialGuess_Update(eps);
340:     }
341:   }
342:   if (!power->update) {
343:     EPSGetStartVector(eps,0,NULL);
344:   }
345:   if (power->nonlinear) {
346:     BVGetColumn(eps->V,0,&v);
347:     EPSPowerUpdateFunctionB(eps,v,Bx);
348:     VecNorm(Bx,NORM_2,&norm);
349:     Normalize(Bx,norm,NULL);
350:     BVRestoreColumn(eps->V,0,&v);
351:   }

353:   STGetShift(eps->st,&sigma);    /* original shift */
354:   rho = sigma;

356:   while (eps->reason == EPS_CONVERGED_ITERATING) {
357:     eps->its++;
358:     k = eps->nconv;

360:     /* y = OP v */
361:     BVGetColumn(eps->V,k,&v);
362:     if (power->nonlinear) {
363:       VecCopy(v,eps->work[3]);
364:       EPSPowerApply_SNES(eps,v,y);
365:       VecCopy(eps->work[3],v);
366:     } else {
367:       STApply(eps->st,v,y);
368:     }
369:     BVRestoreColumn(eps->V,k,&v);

371:     /* purge previously converged eigenvectors */
372:     if (power->nonlinear) {
373:       EPSPowerUpdateFunctionB(eps,y,Bx);
374:       VecNorm(Bx,NORM_2,&norm);
375:       Normalize(Bx,norm,&sign);
376:     } else {
377:       DSGetArray(eps->ds,DS_MAT_A,&T);
378:       BVSetActiveColumns(eps->V,0,k);
379:       BVOrthogonalizeVec(eps->V,y,T+k*ld,&norm,NULL);
380:     }

382:     /* theta = (v,y)_B */
383:     BVSetActiveColumns(eps->V,k,k+1);
384:     BVDotVec(eps->V,y,&theta);
385:     if (!power->nonlinear) {
386:       T[k+k*ld] = theta;
387:       DSRestoreArray(eps->ds,DS_MAT_A,&T);
388:     }

390:     if (power->nonlinear) theta = norm*sign;

392:     if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */

394:       /* approximate eigenvalue is the Rayleigh quotient */
395:       eps->eigr[eps->nconv] = theta;

397:       /* compute relative error as ||y-theta v||_2/|theta| */
398:       VecCopy(y,e);
399:       BVGetColumn(eps->V,k,&v);
400:       VecAXPY(e,power->nonlinear?-1.0:-theta,v);
401:       BVRestoreColumn(eps->V,k,&v);
402:       VecNorm(e,NORM_2,&relerr);
403:       relerr /= PetscAbsScalar(theta);

405:     } else {  /* RQI */

407:       /* delta = ||y||_B */
408:       delta = norm;

410:       /* compute relative error */
411:       if (rho == 0.0) relerr = PETSC_MAX_REAL;
412:       else relerr = 1.0 / (norm*PetscAbsScalar(rho));

414:       /* approximate eigenvalue is the shift */
415:       eps->eigr[eps->nconv] = rho;

417:       /* compute new shift */
418:       if (relerr<eps->tol) {
419:         rho = sigma;  /* if converged, restore original shift */
420:         STSetShift(eps->st,rho);
421:       } else {
422:         rho = rho + theta/(delta*delta);  /* Rayleigh quotient R(v) */
423:         if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
424:           /* beta1 is the norm of the residual associated with R(v) */
425:           BVGetColumn(eps->V,k,&v);
426:           VecAXPY(v,-theta/(delta*delta),y);
427:           BVRestoreColumn(eps->V,k,&v);
428:           BVScaleColumn(eps->V,k,1.0/delta);
429:           BVNormColumn(eps->V,k,NORM_2,&norm1);
430:           beta1 = norm1;

432:           /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
433:           STGetMatrix(eps->st,0,&A);
434:           BVGetColumn(eps->V,k,&v);
435:           MatMult(A,v,e);
436:           VecDot(v,e,&alpha2);
437:           BVRestoreColumn(eps->V,k,&v);
438:           alpha2 = alpha2 / (beta1 * beta1);

440:           /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
441:           PetscFPTrapPush(PETSC_FP_TRAP_OFF);
442:           PetscStackCallBLAS("LAPACKlaev2",LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1));
443:           PetscFPTrapPop();
444:           if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
445:           else rho = rt2;
446:         }
447:         /* update operator according to new shift */
448:         KSPSetErrorIfNotConverged(ksp,PETSC_FALSE);
449:         STSetShift(eps->st,rho);
450:         KSPGetConvergedReason(ksp,&reason);
451:         if (reason) {
452:           PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n");
453:           rho *= 1+PETSC_MACHINE_EPSILON;
454:           STSetShift(eps->st,rho);
455:           KSPGetConvergedReason(ksp,&reason);
456:           if (reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Second factorization failed");
457:         }
458:         KSPSetErrorIfNotConverged(ksp,PETSC_TRUE);
459:       }
460:     }
461:     eps->errest[eps->nconv] = relerr;

463:     /* normalize */
464:     if (!power->nonlinear) { Normalize(y,norm,NULL); }
465:     BVInsertVec(eps->V,k,y);

467:     if (power->update) {
468:       SNESGetConvergedReason(power->snes,&snesreason);
469:     }
470:     /* if relerr<tol, accept eigenpair */
471:     if (relerr<eps->tol || (power->update && snesreason > 0)) {
472:       eps->nconv = eps->nconv + 1;
473:       if (eps->nconv<eps->nev) {
474:         EPSGetStartVector(eps,eps->nconv,&breakdown);
475:         if (breakdown) {
476:           eps->reason = EPS_DIVERGED_BREAKDOWN;
477:           PetscInfo(eps,"Unable to generate more start vectors\n");
478:           break;
479:         }
480:       }
481:     }
482:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev));
483:     (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
484:   }

486:   if (power->nonlinear) {
487:     PetscObjectCompose((PetscObject)power->snes,"eps",NULL);
488:   } else {
489:     DSSetDimensions(eps->ds,eps->nconv,0,0,0);
490:     DSSetState(eps->ds,DS_STATE_RAW);
491:   }
492:   return(0);
493: #endif
494: }


497: PetscErrorCode EPSStopping_Power(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
498: {
500:   EPS_POWER      *power = (EPS_POWER*)eps->data;
501:   SNESConvergedReason snesreason;

504:   if (power->update) {
505:     SNESGetConvergedReason(power->snes,&snesreason);
506:     if (snesreason < 0) {
507:       *reason = EPS_DIVERGED_BREAKDOWN;
508:       return(0);
509:     }
510:   }
511:   EPSStoppingBasic(eps,its,max_it,nconv,nev,reason,ctx);
512:   return(0);
513: }

515: PetscErrorCode EPSBackTransform_Power(EPS eps)
516: {
518:   EPS_POWER      *power = (EPS_POWER*)eps->data;

521:   if (power->nonlinear) eps->eigr[0] = 1.0/eps->eigr[0];
522:   else if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) {
523:     EPSBackTransform_Default(eps);
524:   }
525:   return(0);
526: }

528: PetscErrorCode EPSSetFromOptions_Power(PetscOptionItems *PetscOptionsObject,EPS eps)
529: {
530:   PetscErrorCode    ierr;
531:   EPS_POWER         *power = (EPS_POWER*)eps->data;
532:   PetscBool         flg,val;
533:   EPSPowerShiftType shift;

536:   PetscOptionsHead(PetscOptionsObject,"EPS Power Options");

538:     PetscOptionsEnum("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",EPSPowerShiftTypes,(PetscEnum)power->shift_type,(PetscEnum*)&shift,&flg);
539:     if (flg) { EPSPowerSetShiftType(eps,shift); }

541:     PetscOptionsBool("-eps_power_nonlinear","Use nonlinear inverse iteration","EPSPowerSetNonlinear",power->nonlinear,&val,&flg);
542:     if (flg) { EPSPowerSetNonlinear(eps,val); }

544:     PetscOptionsBool("-eps_power_update","Update residual monolithically","EPSPowerSetUpdate",power->update,&val,&flg);
545:     if (flg) { EPSPowerSetUpdate(eps,val); }

547:   PetscOptionsTail();
548:   return(0);
549: }

551: static PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
552: {
553:   EPS_POWER *power = (EPS_POWER*)eps->data;

556:   switch (shift) {
557:     case EPS_POWER_SHIFT_CONSTANT:
558:     case EPS_POWER_SHIFT_RAYLEIGH:
559:     case EPS_POWER_SHIFT_WILKINSON:
560:       if (power->shift_type != shift) {
561:         power->shift_type = shift;
562:         eps->state = EPS_STATE_INITIAL;
563:       }
564:       break;
565:     default:
566:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
567:   }
568:   return(0);
569: }

571: /*@
572:    EPSPowerSetShiftType - Sets the type of shifts used during the power
573:    iteration. This can be used to emulate the Rayleigh Quotient Iteration
574:    (RQI) method.

576:    Logically Collective on EPS

578:    Input Parameters:
579: +  eps - the eigenproblem solver context
580: -  shift - the type of shift

582:    Options Database Key:
583: .  -eps_power_shift_type - Sets the shift type (either 'constant' or
584:                            'rayleigh' or 'wilkinson')

586:    Notes:
587:    By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
588:    is the simple power method (or inverse iteration if a shift-and-invert
589:    transformation is being used).

591:    A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
592:    EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
593:    a cubic converging method such as RQI.

595:    Level: advanced

597: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
598: @*/
599: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
600: {

606:   PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
607:   return(0);
608: }

610: static PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
611: {
612:   EPS_POWER *power = (EPS_POWER*)eps->data;

615:   *shift = power->shift_type;
616:   return(0);
617: }

619: /*@
620:    EPSPowerGetShiftType - Gets the type of shifts used during the power
621:    iteration.

623:    Not Collective

625:    Input Parameter:
626: .  eps - the eigenproblem solver context

628:    Input Parameter:
629: .  shift - the type of shift

631:    Level: advanced

633: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
634: @*/
635: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
636: {

642:   PetscUseMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
643:   return(0);
644: }

646: static PetscErrorCode EPSPowerSetNonlinear_Power(EPS eps,PetscBool nonlinear)
647: {
648:   EPS_POWER *power = (EPS_POWER*)eps->data;

651:   if (power->nonlinear != nonlinear) {
652:     power->nonlinear = nonlinear;
653:     eps->useds = PetscNot(nonlinear);
654:     eps->state = EPS_STATE_INITIAL;
655:   }
656:   return(0);
657: }

659: /*@
660:    EPSPowerSetNonlinear - Sets a flag to indicate that the problem is nonlinear.

662:    Logically Collective on EPS

664:    Input Parameters:
665: +  eps - the eigenproblem solver context
666: -  nonlinear - whether the problem is nonlinear or not

668:    Options Database Key:
669: .  -eps_power_nonlinear - Sets the nonlinear flag

671:    Notes:
672:    If this flag is set, the solver assumes that the problem is nonlinear,
673:    that is, the operators that define the eigenproblem are not constant
674:    matrices, but depend on the eigenvector: A(x)*x=lambda*B(x)*x. This is
675:    different from the case of nonlinearity with respect to the eigenvalue
676:    (use the NEP solver class for this kind of problems).

678:    The way in which nonlinear operators are specified is very similar to
679:    the case of PETSc's SNES solver. The difference is that the callback
680:    functions are provided via composed functions "formFunction" and
681:    "formJacobian" in each of the matrix objects passed as arguments of
682:    EPSSetOperators(). The application context required for these functions
683:    can be attached via a composed PetscContainer.

685:    Level: advanced

687: .seealso: EPSPowerGetNonlinear(), EPSSetOperators()
688: @*/
689: PetscErrorCode EPSPowerSetNonlinear(EPS eps,PetscBool nonlinear)
690: {

696:   PetscTryMethod(eps,"EPSPowerSetNonlinear_C",(EPS,PetscBool),(eps,nonlinear));
697:   return(0);
698: }

700: static PetscErrorCode EPSPowerGetNonlinear_Power(EPS eps,PetscBool *nonlinear)
701: {
702:   EPS_POWER *power = (EPS_POWER*)eps->data;

705:   *nonlinear = power->nonlinear;
706:   return(0);
707: }

709: /*@
710:    EPSPowerGetNonlinear - Returns a flag indicating if the problem is nonlinear.

712:    Not Collective

714:    Input Parameter:
715: .  eps - the eigenproblem solver context

717:    Input Parameter:
718: .  nonlinear - the nonlinear flag

720:    Level: advanced

722: .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
723: @*/
724: PetscErrorCode EPSPowerGetNonlinear(EPS eps,PetscBool *nonlinear)
725: {

731:   PetscUseMethod(eps,"EPSPowerGetNonlinear_C",(EPS,PetscBool*),(eps,nonlinear));
732:   return(0);
733: }

735: static PetscErrorCode EPSPowerSetUpdate_Power(EPS eps,PetscBool update)
736: {
737:   EPS_POWER *power = (EPS_POWER*)eps->data;

740:   if (!power->nonlinear) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"This option does not make sense for linear problems \n");
741:   power->update = update;
742:   eps->state = EPS_STATE_INITIAL;
743:   return(0);
744: }

746: /*@
747:    EPSPowerSetUpdate - Sets a flag to indicate that the residual is updated monolithically
748:    for nonlinear problems. This potentially has a better convergence.

750:    Logically Collective on EPS

752:    Input Parameters:
753: +  eps - the eigenproblem solver context
754: -  update - whether the residual is updated monolithically or not

756:    Options Database Key:
757: .  -eps_power_update - Sets the update flag

759:    Level: advanced

761: .seealso: EPSPowerGetUpdate(), EPSPowerGetNonlinear(), EPSSetOperators()
762: @*/
763: PetscErrorCode EPSPowerSetUpdate(EPS eps,PetscBool update)
764: {

770:   PetscTryMethod(eps,"EPSPowerSetUpdate_C",(EPS,PetscBool),(eps,update));
771:   return(0);
772: }

774: static PetscErrorCode EPSPowerGetUpdate_Power(EPS eps,PetscBool *update)
775: {
776:   EPS_POWER *power = (EPS_POWER*)eps->data;

779:   *update = power->update;
780:   return(0);
781: }

783: /*@
784:    EPSPowerGetUpdate - Returns a flag indicating if the residual is updated monolithically
785:    for nonlinear problems.

787:    Not Collective

789:    Input Parameter:
790: .  eps - the eigenproblem solver context

792:    Input Parameter:
793: .  update - the update flag

795:    Level: advanced

797: .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
798: @*/
799: PetscErrorCode EPSPowerGetUpdate(EPS eps,PetscBool *update)
800: {

806:   PetscUseMethod(eps,"EPSPowerGetUpdate_C",(EPS,PetscBool*),(eps,update));
807:   return(0);
808: }

810: static PetscErrorCode EPSPowerSetSNES_Power(EPS eps,SNES snes)
811: {
813:   EPS_POWER      *power = (EPS_POWER*)eps->data;

816:   PetscObjectReference((PetscObject)snes);
817:   SNESDestroy(&power->snes);
818:   power->snes = snes;
819:   PetscLogObjectParent((PetscObject)eps,(PetscObject)power->snes);
820:   eps->state = EPS_STATE_INITIAL;
821:   return(0);
822: }

824: /*@
825:    EPSPowerSetSNES - Associate a nonlinear solver object (SNES) to the
826:    eigenvalue solver (to be used in nonlinear inverse iteration).

828:    Collective on EPS

830:    Input Parameters:
831: +  eps  - the eigenvalue solver
832: -  snes - the nonlinear solver object

834:    Level: advanced

836: .seealso: EPSPowerGetSNES()
837: @*/
838: PetscErrorCode EPSPowerSetSNES(EPS eps,SNES snes)
839: {

846:   PetscTryMethod(eps,"EPSPowerSetSNES_C",(EPS,SNES),(eps,snes));
847:   return(0);
848: }

850: static PetscErrorCode EPSPowerGetSNES_Power(EPS eps,SNES *snes)
851: {
853:   EPS_POWER      *power = (EPS_POWER*)eps->data;

856:   if (!power->snes) {
857:     SNESCreate(PetscObjectComm((PetscObject)eps),&power->snes);
858:     SNESSetOptionsPrefix(power->snes,((PetscObject)eps)->prefix);
859:     SNESAppendOptionsPrefix(power->snes,"eps_power_");
860:     PetscObjectIncrementTabLevel((PetscObject)power->snes,(PetscObject)eps,1);
861:     PetscLogObjectParent((PetscObject)eps,(PetscObject)power->snes);
862:   }
863:   *snes = power->snes;
864:   return(0);
865: }

867: /*@
868:    EPSPowerGetSNES - Retrieve the nonlinear solver object (SNES) associated
869:    with the eigenvalue solver.

871:    Not Collective

873:    Input Parameter:
874: .  eps - the eigenvalue solver

876:    Output Parameter:
877: .  snes - the nonlinear solver object

879:    Level: advanced

881: .seealso: EPSPowerSetSNES()
882: @*/
883: PetscErrorCode EPSPowerGetSNES(EPS eps,SNES *snes)
884: {

890:   PetscUseMethod(eps,"EPSPowerGetSNES_C",(EPS,SNES*),(eps,snes));
891:   return(0);
892: }

894: PetscErrorCode EPSReset_Power(EPS eps)
895: {
897:   EPS_POWER      *power = (EPS_POWER*)eps->data;

900:   if (power->snes) { SNESReset(power->snes); }
901:   return(0);
902: }

904: PetscErrorCode EPSDestroy_Power(EPS eps)
905: {
907:   EPS_POWER      *power = (EPS_POWER*)eps->data;

910:   if (power->nonlinear) {
911:     SNESDestroy(&power->snes);
912:   }
913:   PetscFree(eps->data);
914:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",NULL);
915:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",NULL);
916:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",NULL);
917:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",NULL);
918:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",NULL);
919:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",NULL);
920:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",NULL);
921:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",NULL);
922:   return(0);
923: }

925: PetscErrorCode EPSView_Power(EPS eps,PetscViewer viewer)
926: {
928:   EPS_POWER      *power = (EPS_POWER*)eps->data;
929:   PetscBool      isascii;

932:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
933:   if (isascii) {
934:     if (power->nonlinear) {
935:       PetscViewerASCIIPrintf(viewer,"  using nonlinear inverse iteration\n");
936:       if (power->update) {
937:         PetscViewerASCIIPrintf(viewer,"  updating the residual monolithically\n");
938:       }
939:       if (!power->snes) { EPSPowerGetSNES(eps,&power->snes); }
940:       PetscViewerASCIIPushTab(viewer);
941:       SNESView(power->snes,viewer);
942:       PetscViewerASCIIPopTab(viewer);
943:     } else {
944:       PetscViewerASCIIPrintf(viewer,"  %s shifts\n",EPSPowerShiftTypes[power->shift_type]);
945:     }
946:   }
947:   return(0);
948: }

950: PetscErrorCode EPSComputeVectors_Power(EPS eps)
951: {
953:   EPS_POWER      *power = (EPS_POWER*)eps->data;

956:   if (!power->nonlinear) {
957:     EPSComputeVectors_Schur(eps);
958:   }
959:   return(0);
960: }

962: PetscErrorCode EPSSetDefaultST_Power(EPS eps)
963: {
965:   EPS_POWER      *power = (EPS_POWER*)eps->data;
966:   KSP            ksp;
967:   PC             pc;

970:   if (power->nonlinear) {
971:     STGetKSP(eps->st,&ksp);
972:     KSPGetPC(ksp,&pc);
973:     PCSetType(pc,PCNONE);
974:   }
975:   return(0);
976: }

978: PETSC_EXTERN PetscErrorCode EPSCreate_Power(EPS eps)
979: {
980:   EPS_POWER      *ctx;

984:   PetscNewLog(eps,&ctx);
985:   eps->data = (void*)ctx;

987:   eps->useds = PETSC_TRUE;
988:   eps->categ = EPS_CATEGORY_OTHER;

990:   eps->ops->solve          = EPSSolve_Power;
991:   eps->ops->setup          = EPSSetUp_Power;
992:   eps->ops->setfromoptions = EPSSetFromOptions_Power;
993:   eps->ops->reset          = EPSReset_Power;
994:   eps->ops->destroy        = EPSDestroy_Power;
995:   eps->ops->view           = EPSView_Power;
996:   eps->ops->backtransform  = EPSBackTransform_Power;
997:   eps->ops->computevectors = EPSComputeVectors_Power;
998:   eps->ops->setdefaultst   = EPSSetDefaultST_Power;
999:   eps->stopping            = EPSStopping_Power;

1001:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",EPSPowerSetShiftType_Power);
1002:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",EPSPowerGetShiftType_Power);
1003:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",EPSPowerSetNonlinear_Power);
1004:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",EPSPowerGetNonlinear_Power);
1005:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",EPSPowerSetUpdate_Power);
1006:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",EPSPowerGetUpdate_Power);
1007:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",EPSPowerSetSNES_Power);
1008:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",EPSPowerGetSNES_Power);
1009:   return(0);
1010: }