Actual source code: pepsolve.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:       PEP routines related to the solution process.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

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

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

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

 24: #include <slepc-private/pepimpl.h>       /*I "slepcpep.h" I*/
 25: #include <petscdraw.h>

 29: PETSC_STATIC_INLINE PetscErrorCode PEPComputeVectors(PEP pep)
 30: {

 34:   PEPCheckSolved(pep,1);
 35:   switch (pep->state) {
 36:   case PEP_STATE_SOLVED:
 37:     if (pep->ops->computevectors) {
 38:       (*pep->ops->computevectors)(pep);
 39:     }
 40:     break;
 41:   default:
 42:     break;
 43:   }
 44:   pep->state = PEP_STATE_EIGENVECTORS;
 45:   return(0);
 46: }

 50: /*@
 51:    PEPSolve - Solves the polynomial eigensystem.

 53:    Collective on PEP

 55:    Input Parameter:
 56: .  pep - eigensolver context obtained from PEPCreate()

 58:    Options Database Keys:
 59: +  -pep_view - print information about the solver used
 60: -  -pep_plot_eigs - plot computed eigenvalues

 62:    Level: beginner

 64: .seealso: PEPCreate(), PEPSetUp(), PEPDestroy(), PEPSetTolerances()
 65: @*/
 66: PetscErrorCode PEPSolve(PEP pep)
 67: {
 68:   PetscErrorCode    ierr;
 69:   PetscInt          i;
 70:   PetscReal         re,im;
 71:   PetscBool         flg,islinear;
 72:   PetscViewer       viewer;
 73:   PetscViewerFormat format;
 74:   PetscDraw         draw;
 75:   PetscDrawSP       drawsp;

 79:   PetscLogEventBegin(PEP_Solve,pep,0,0,0);

 81:   /* call setup */
 82:   PEPSetUp(pep);
 83:   pep->nconv = 0;
 84:   pep->its   = 0;
 85:   for (i=0;i<pep->ncv;i++) {
 86:     pep->eigr[i]   = 0.0;
 87:     pep->eigi[i]   = 0.0;
 88:     pep->errest[i] = 0.0;
 89:   }
 90:   PEPMonitor(pep,pep->its,pep->nconv,pep->eigr,pep->eigi,pep->errest,pep->ncv);

 92:   (*pep->ops->solve)(pep);
 93:   
 94:   PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&islinear);
 95:   if (!islinear) {
 96:     STPostSolve(pep->st);
 97:   }

 99:   if (!pep->reason) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");

101:   if (!islinear) {
102:     /* Map eigenvalues back to the original problem */
103:     STGetTransform(pep->st,&flg);
104:     if (flg) {
105:       STBackTransform(pep->st,pep->nconv,pep->eigr,pep->eigi);
106:     }
107:   }

109:   pep->state = PEP_STATE_SOLVED;

111:   if (pep->refine==PEP_REFINE_SIMPLE && pep->rits>0) {
112:     PEPComputeVectors(pep);
113:     PEPNewtonRefinementSimple(pep,&pep->rits,&pep->rtol,pep->nconv);
114:     pep->state = PEP_STATE_EIGENVECTORS;
115:   }

117: #if !defined(PETSC_USE_COMPLEX)
118:   /* reorder conjugate eigenvalues (positive imaginary first) */
119:   for (i=0;i<pep->nconv-1;i++) {
120:     if (pep->eigi[i] != 0) {
121:       if (pep->eigi[i] < 0) {
122:         pep->eigi[i] = -pep->eigi[i];
123:         pep->eigi[i+1] = -pep->eigi[i+1];
124:         /* the next correction only works with eigenvectors */
125:         PEPComputeVectors(pep);
126:         BVScaleColumn(pep->V,i+1,-1.0);
127:       }
128:       i++;
129:     }
130:   }
131: #endif

133:   /* sort eigenvalues according to pep->which parameter */
134:   SlepcSortEigenvalues(pep->sc,pep->nconv,pep->eigr,pep->eigi,pep->perm);

136:   PetscLogEventEnd(PEP_Solve,pep,0,0,0);

138:   /* various viewers */
139:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,"-pep_view",&viewer,&format,&flg);
140:   if (flg && !PetscPreLoadingOn) {
141:     PetscViewerPushFormat(viewer,format);
142:     PEPView(pep,viewer);
143:     PetscViewerPopFormat(viewer);
144:     PetscViewerDestroy(&viewer);
145:   }

147:   flg = PETSC_FALSE;
148:   PetscOptionsGetBool(((PetscObject)pep)->prefix,"-pep_plot_eigs",&flg,NULL);
149:   if (flg) {
150:     PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
151:     PetscViewerDrawGetDraw(viewer,0,&draw);
152:     PetscDrawSPCreate(draw,1,&drawsp);
153:     for (i=0;i<pep->nconv;i++) {
154: #if defined(PETSC_USE_COMPLEX)
155:       re = PetscRealPart(pep->eigr[i]);
156:       im = PetscImaginaryPart(pep->eigi[i]);
157: #else
158:       re = pep->eigr[i];
159:       im = pep->eigi[i];
160: #endif
161:       PetscDrawSPAddPoint(drawsp,&re,&im);
162:     }
163:     PetscDrawSPDraw(drawsp,PETSC_TRUE);
164:     PetscDrawSPDestroy(&drawsp);
165:     PetscViewerDestroy(&viewer);
166:   }

168:   /* Remove the initial subspace */
169:   pep->nini = 0;
170:   return(0);
171: }

175: /*@
176:    PEPGetIterationNumber - Gets the current iteration number. If the
177:    call to PEPSolve() is complete, then it returns the number of iterations
178:    carried out by the solution method.

180:    Not Collective

182:    Input Parameter:
183: .  pep - the polynomial eigensolver context

185:    Output Parameter:
186: .  its - number of iterations

188:    Level: intermediate

190:    Note:
191:    During the i-th iteration this call returns i-1. If PEPSolve() is
192:    complete, then parameter "its" contains either the iteration number at
193:    which convergence was successfully reached, or failure was detected.
194:    Call PEPGetConvergedReason() to determine if the solver converged or
195:    failed and why.

197: .seealso: PEPGetConvergedReason(), PEPSetTolerances()
198: @*/
199: PetscErrorCode PEPGetIterationNumber(PEP pep,PetscInt *its)
200: {
204:   *its = pep->its;
205:   return(0);
206: }

210: /*@
211:    PEPGetConverged - Gets the number of converged eigenpairs.

213:    Not Collective

215:    Input Parameter:
216: .  pep - the polynomial eigensolver context

218:    Output Parameter:
219: .  nconv - number of converged eigenpairs

221:    Note:
222:    This function should be called after PEPSolve() has finished.

224:    Level: beginner

226: .seealso: PEPSetDimensions(), PEPSolve()
227: @*/
228: PetscErrorCode PEPGetConverged(PEP pep,PetscInt *nconv)
229: {
233:   PEPCheckSolved(pep,1);
234:   *nconv = pep->nconv;
235:   return(0);
236: }

240: /*@C
241:    PEPGetConvergedReason - Gets the reason why the PEPSolve() iteration was
242:    stopped.

244:    Not Collective

246:    Input Parameter:
247: .  pep - the polynomial eigensolver context

249:    Output Parameter:
250: .  reason - negative value indicates diverged, positive value converged

252:    Possible values for reason:
253: +  PEP_CONVERGED_TOL - converged up to tolerance
254: .  PEP_DIVERGED_ITS - required more than its to reach convergence
255: -  PEP_DIVERGED_BREAKDOWN - generic breakdown in method

257:    Note:
258:    Can only be called after the call to PEPSolve() is complete.

260:    Level: intermediate

262: .seealso: PEPSetTolerances(), PEPSolve(), PEPConvergedReason
263: @*/
264: PetscErrorCode PEPGetConvergedReason(PEP pep,PEPConvergedReason *reason)
265: {
269:   PEPCheckSolved(pep,1);
270:   *reason = pep->reason;
271:   return(0);
272: }

276: /*@
277:    PEPGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
278:    PEPSolve(). The solution consists in both the eigenvalue and the eigenvector.

280:    Logically Collective on EPS

282:    Input Parameters:
283: +  pep - polynomial eigensolver context
284: -  i   - index of the solution

286:    Output Parameters:
287: +  eigr - real part of eigenvalue
288: .  eigi - imaginary part of eigenvalue
289: .  Vr   - real part of eigenvector
290: -  Vi   - imaginary part of eigenvector

292:    Notes:
293:    If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
294:    configured with complex scalars the eigenvalue is stored
295:    directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
296:    set to zero).

298:    The index i should be a value between 0 and nconv-1 (see PEPGetConverged()).
299:    Eigenpairs are indexed according to the ordering criterion established
300:    with PEPSetWhichEigenpairs().

302:    Level: beginner

304: .seealso: PEPSolve(), PEPGetConverged(), PEPSetWhichEigenpairs()
305: @*/
306: PetscErrorCode PEPGetEigenpair(PEP pep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
307: {
308:   PetscInt       k;

316:   PEPCheckSolved(pep,1);
317:   if (i<0 || i>=pep->nconv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");

319:   PEPComputeVectors(pep);
320:   if (!pep->perm) k = i;
321:   else k = pep->perm[i];

323:   /* eigenvalue */
324: #if defined(PETSC_USE_COMPLEX)
325:   if (eigr) *eigr = pep->eigr[k];
326:   if (eigi) *eigi = 0;
327: #else
328:   if (eigr) *eigr = pep->eigr[k];
329:   if (eigi) *eigi = pep->eigi[k];
330: #endif

332:   /* eigenvector */
333: #if defined(PETSC_USE_COMPLEX)
334:   if (Vr) { BVCopyVec(pep->V,k,Vr); }
335:   if (Vi) { VecSet(Vi,0.0); }
336: #else
337:   if (pep->eigi[k]>0) { /* first value of conjugate pair */
338:     if (Vr) { BVCopyVec(pep->V,k,Vr); }
339:     if (Vi) { BVCopyVec(pep->V,k+1,Vi); }
340:   } else if (pep->eigi[k]<0) { /* second value of conjugate pair */
341:     if (Vr) { BVCopyVec(pep->V,k-1,Vr); }
342:     if (Vi) {
343:       BVCopyVec(pep->V,k,Vi);
344:       VecScale(Vi,-1.0);
345:     }
346:   } else { /* real eigenvalue */
347:     if (Vr) { BVCopyVec(pep->V,k,Vr); }
348:     if (Vi) { VecSet(Vi,0.0); }
349:   }
350: #endif
351:   return(0);
352: }

356: /*@
357:    PEPGetErrorEstimate - Returns the error estimate associated to the i-th
358:    computed eigenpair.

360:    Not Collective

362:    Input Parameter:
363: +  pep - polynomial eigensolver context
364: -  i   - index of eigenpair

366:    Output Parameter:
367: .  errest - the error estimate

369:    Notes:
370:    This is the error estimate used internally by the eigensolver. The actual
371:    error bound can be computed with PEPComputeRelativeError(). See also the users
372:    manual for details.

374:    Level: advanced

376: .seealso: PEPComputeRelativeError()
377: @*/
378: PetscErrorCode PEPGetErrorEstimate(PEP pep,PetscInt i,PetscReal *errest)
379: {
383:   PEPCheckSolved(pep,1);
384:   if (i<0 || i>=pep->nconv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
385:   if (pep->perm) i = pep->perm[i];
386:   if (errest) *errest = pep->errest[i];
387:   return(0);
388: }

392: /*
393:    PEPComputeResidualNorm_Private - Computes the norm of the residual vector
394:    associated with an eigenpair.
395: */
396: PetscErrorCode PEPComputeResidualNorm_Private(PEP pep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *norm)
397: {
399:   Vec            u,w;
400:   Mat            *A=pep->A;
401:   PetscInt       i,nmat=pep->nmat;
402:   PetscScalar    t[20],*vals=t,*ivals=NULL;
403: #if !defined(PETSC_USE_COMPLEX)
404:   Vec            ui,wi;
405:   PetscReal      ni;
406:   PetscBool      imag;
407:   PetscScalar    it[20];
408: #endif

411:   BVGetVec(pep->V,&u);
412:   BVGetVec(pep->V,&w);
413:   VecZeroEntries(u);
414: #if !defined(PETSC_USE_COMPLEX)
415:   ivals = it; 
416: #endif
417:   if (nmat>20) {
418:     PetscMalloc(nmat*sizeof(PetscScalar),&vals);
419: #if !defined(PETSC_USE_COMPLEX)
420:     PetscMalloc(nmat*sizeof(PetscScalar),&ivals);
421: #endif
422:   }
423:   PEPEvaluateBasis(pep,kr,ki,vals,ivals);
424: #if !defined(PETSC_USE_COMPLEX)
425:   if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON))
426:     imag = PETSC_FALSE;
427:   else {
428:     imag = PETSC_TRUE;
429:     VecDuplicate(u,&ui);
430:     VecDuplicate(u,&wi);
431:     VecZeroEntries(ui);
432:   }
433: #endif
434:   for (i=0;i<nmat;i++) {
435:     if (vals[i]!=0.0) {
436:       MatMult(A[i],xr,w);
437:       VecAXPY(u,vals[i],w);
438:     }
439: #if !defined(PETSC_USE_COMPLEX)
440:     if (imag) {
441:       if (ivals[i]!=0 || vals[i]!=0) {
442:         MatMult(A[i],xi,wi);
443:         if (vals[i]==0) {
444:           MatMult(A[i],xr,w);
445:         }
446:       }
447:       if (ivals[i]!=0){
448:         VecAXPY(u,-ivals[i],wi);
449:         VecAXPY(ui,ivals[i],w);
450:       }
451:       if (vals[i]!=0) {
452:         VecAXPY(ui,vals[i],wi);
453:       }
454:     }
455: #endif
456:   }
457:   VecNorm(u,NORM_2,norm);
458: #if !defined(PETSC_USE_COMPLEX)
459:   if (imag) {
460:     VecNorm(ui,NORM_2,&ni);
461:     *norm = SlepcAbsEigenvalue(*norm,ni);
462:   }
463: #endif
464:   VecDestroy(&w);
465:   VecDestroy(&u);
466:   if (nmat>20) {
467:     PetscFree(vals);
468: #if !defined(PETSC_USE_COMPLEX)
469:     PetscFree(ivals);
470: #endif
471:   }
472: #if !defined(PETSC_USE_COMPLEX)
473:   if (imag) {
474:     VecDestroy(&wi);
475:     VecDestroy(&ui);
476:   }
477: #endif
478:   return(0);
479: }

483: /*@
484:    PEPComputeResidualNorm - Computes the norm of the residual vector associated
485:    with the i-th computed eigenpair.

487:    Collective on PEP

489:    Input Parameter:
490: +  pep - the polynomial eigensolver context
491: -  i   - the solution index

493:    Output Parameter:
494: .  norm - the residual norm, computed as ||P(l)x||_2 where l is the
495:    eigenvalue and x is the eigenvector.

497:    Notes:
498:    The index i should be a value between 0 and nconv-1 (see PEPGetConverged()).
499:    Eigenpairs are indexed according to the ordering criterion established
500:    with PEPSetWhichEigenpairs().

502:    Level: beginner

504: .seealso: PEPSolve(), PEPGetConverged(), PEPSetWhichEigenpairs()
505: @*/
506: PetscErrorCode PEPComputeResidualNorm(PEP pep,PetscInt i,PetscReal *norm)
507: {
509:   Vec            xr,xi;
510:   PetscScalar    kr,ki;

516:   PEPCheckSolved(pep,1);
517:   BVGetVec(pep->V,&xr);
518:   BVGetVec(pep->V,&xi);
519:   PEPGetEigenpair(pep,i,&kr,&ki,xr,xi);
520:   PEPComputeResidualNorm_Private(pep,kr,ki,xr,xi,norm);
521:   VecDestroy(&xr);
522:   VecDestroy(&xi);
523:   return(0);
524: }

528: /*
529:    PEPComputeRelativeError_Private - Computes the relative error bound
530:    associated with an eigenpair.
531: */
532: PetscErrorCode PEPComputeRelativeError_Private(PEP pep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *error)
533: {
535:   PetscReal      norm,er;
536: #if !defined(PETSC_USE_COMPLEX)
537:   PetscReal      ei;
538: #endif

541:   PEPComputeResidualNorm_Private(pep,kr,ki,xr,xi,&norm);
542: #if !defined(PETSC_USE_COMPLEX)
543:   if (ki == 0) {
544: #endif
545:     VecNorm(xr,NORM_2,&er);
546: #if !defined(PETSC_USE_COMPLEX)
547:   } else {
548:     VecNorm(xr,NORM_2,&er);
549:     VecNorm(xi,NORM_2,&ei);
550:     er = SlepcAbsEigenvalue(er,ei);
551:   }
552: #endif
553:   (*pep->converged)(pep,kr,ki,norm/er,error,pep->convergedctx);
554:   return(0);
555: }

559: /*@
560:    PEPComputeRelativeError - Computes the relative error bound associated
561:    with the i-th computed eigenpair.

563:    Collective on PEP

565:    Input Parameter:
566: +  pep - the polynomial eigensolver context
567: -  i   - the solution index

569:    Output Parameter:
570: .  error - the relative error bound, computed as ||P(l)x||_2/||lx||_2 where
571:    l is the eigenvalue and x is the eigenvector.

573:    Level: beginner

575: .seealso: PEPSolve(), PEPComputeResidualNorm(), PEPGetErrorEstimate()
576: @*/
577: PetscErrorCode PEPComputeRelativeError(PEP pep,PetscInt i,PetscReal *error)
578: {
580:   Vec            xr,xi;
581:   PetscScalar    kr,ki;

587:   PEPCheckSolved(pep,1);
588:   BVGetVec(pep->V,&xr);
589:   BVGetVec(pep->V,&xi);
590:   PEPGetEigenpair(pep,i,&kr,&ki,xr,xi);
591:   PEPComputeRelativeError_Private(pep,kr,ki,xr,xi,error);
592:   VecDestroy(&xr);
593:   VecDestroy(&xi);
594:   return(0);
595: }