Actual source code: solve.c

  1: /*
  2:       EPS routines related to the solution process.
  3: */
 4:  #include src/eps/epsimpl.h

  8: /*@
  9:    EPSSolve - Solves the eigensystem.

 11:    Collective on EPS

 13:    Input Parameter:
 14: .  eps - eigensolver context obtained from EPSCreate()

 16:    Options Database:
 17: +   -eps_view - print information about the solver used
 18: .   -eps_view_binary - save the matrices to the default binary file
 19: -   -eps_plot_eigs - plot computed eigenvalues

 21:    Level: beginner

 23: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances() 
 24: @*/
 25: PetscErrorCode EPSSolve(EPS eps)
 26: {
 27:   PetscErrorCode ierr;
 28:   int            i;
 29:   PetscReal      re,im;
 30:   PetscTruth     flg;
 31:   PetscViewer    viewer;
 32:   PetscDraw      draw;
 33:   PetscDrawSP    drawsp;


 38:   PetscOptionsHasName(eps->prefix,"-eps_view_binary",&flg);
 39:   if (flg) {
 40:     Mat A,B;
 41:     STGetOperators(eps->OP,&A,&B);
 42:     MatView(A,PETSC_VIEWER_BINARY_(eps->comm));
 43:     if (B) MatView(B,PETSC_VIEWER_BINARY_(eps->comm));
 44:   }

 46:   /* reset the convergence flag from the previous solves */
 47:   eps->reason = EPS_CONVERGED_ITERATING;

 49:   if (!eps->setupcalled){ EPSSetUp(eps); }
 50:   STResetNumberLinearIterations(eps->OP);
 51:   eps->evecsavailable = PETSC_FALSE;
 52:   PetscLogEventBegin(EPS_Solve,eps,eps->V[0],eps->V[0],0);
 53:   STPreSolve(eps->OP,eps);
 54:   (*eps->ops->solve)(eps);
 55:   STPostSolve(eps->OP,eps);
 56:   PetscLogEventEnd(EPS_Solve,eps,eps->V[0],eps->V[0],0);
 57:   if (!eps->reason) {
 58:     SETERRQ(1,"Internal error, solver returned without setting converged reason");
 59:   }

 61:   /* Map eigenvalues back to the original problem, necessary in some 
 62:   * spectral transformations */
 63:   (*eps->ops->backtransform)(eps);

 65: #ifndef PETSC_USE_COMPLEX
 66:   /* reorder conjugate eigenvalues (positive imaginary first) */
 67:   for (i=0; i<eps->nconv-1; i++) {
 68:     PetscScalar minus = -1.0;
 69:     if (eps->eigi[i] != 0) {
 70:       if (eps->eigi[i] < 0) {
 71:         eps->eigi[i] = -eps->eigi[i];
 72:         eps->eigi[i+1] = -eps->eigi[i+1];
 73:         VecScale(&minus, eps->V[i+1]);
 74:       }
 75:       i++;
 76:     }
 77:   }
 78: #endif

 80:   /* sort eigenvalues according to eps->which parameter */
 81:   if (eps->perm) {
 82:     PetscFree(eps->perm);
 83:     eps->perm = PETSC_NULL;
 84:   }
 85:   if (eps->nconv > 0) {
 86:     PetscMalloc(sizeof(int)*eps->nconv, &eps->perm);
 87:     EPSSortEigenvalues(eps->nconv, eps->eigr, eps->eigi, eps->which, eps->nconv, eps->perm);
 88:   }

 90:   PetscOptionsHasName(eps->prefix,"-eps_view",&flg);
 91:   if (flg && !PetscPreLoadingOn) { EPSView(eps,PETSC_VIEWER_STDOUT_WORLD); }

 93:   PetscOptionsHasName(eps->prefix,"-eps_plot_eigs",&flg);
 94:   if (flg) {
 95:     PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
 96:                              PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
 97:     PetscViewerDrawGetDraw(viewer,0,&draw);
 98:     PetscDrawSPCreate(draw,1,&drawsp);
 99:     for( i=0; i<eps->nconv; i++ ) {
100: #if defined(PETSC_USE_COMPLEX)
101:       re = PetscRealPart(eps->eigr[i]);
102:       im = PetscImaginaryPart(eps->eigi[i]);
103: #else
104:       re = eps->eigr[i];
105:       im = eps->eigi[i];
106: #endif
107:       PetscDrawSPAddPoint(drawsp,&re,&im);
108:     }
109:     PetscDrawSPDraw(drawsp);
110:     PetscDrawSPDestroy(drawsp);
111:     PetscViewerDestroy(viewer);
112:   }

114:   return(0);
115: }

119: /*@
120:    EPSGetIterationNumber - Gets the current iteration number. If the 
121:    call to EPSSolve() is complete, then it returns the number of iterations 
122:    carried out by the solution method.
123:  
124:    Not Collective

126:    Input Parameter:
127: .  eps - the eigensolver context

129:    Output Parameter:
130: .  its - number of iterations

132:    Level: intermediate

134:    Notes:
135:       During the i-th iteration this call returns i-1. If EPSSolve() is 
136:       complete, then parameter "its" contains either the iteration number at
137:       which convergence was successfully reached, or failure was detected.  
138:       Call EPSGetConvergedReason() to determine if the solver converged or 
139:       failed and why.

141: @*/
142: PetscErrorCode EPSGetIterationNumber(EPS eps,int *its)
143: {
147:   *its = eps->its;
148:   return(0);
149: }

153: /*@
154:    EPSGetNumberLinearIterations - Gets the total number of iterations
155:    required by the linear solves associated to the ST object during the 
156:    last EPSSolve() call.

158:    Not Collective

160:    Input Parameter:
161: .  eps - EPS context

163:    Output Parameter:
164: .  lits - number of linear iterations

166:    Notes:
167:    When the eigensolver algorithm invokes STApply() then a linear system 
168:    must be solved (except in the case of standard eigenproblems and shift
169:    transformation). The number of iterations required in this solve is
170:    accumulated into a counter whose value is returned by this function.

172:    The iteration counter is reset to zero at each successive call to EPSSolve().

174:    Level: intermediate

176: @*/
177: PetscErrorCode EPSGetNumberLinearIterations(EPS eps,int* lits)
178: {
182:   STGetNumberLinearIterations(eps->OP, lits);
183:   return(0);
184: }

188: /*@
189:    EPSGetConverged - Gets the number of converged eigenpairs.

191:    Not Collective

193:    Input Parameter:
194: .  eps - the eigensolver context
195:   
196:    Output Parameter:
197: .  nconv - number of converged eigenpairs 

199:    Note:
200:    This function should be called after EPSSolve() has finished.

202:    Level: beginner

204: .seealso: EPSSetDimensions()
205: @*/
206: PetscErrorCode EPSGetConverged(EPS eps,int *nconv)
207: {
210:   if (nconv) *nconv = eps->nconv;
211:   return(0);
212: }


217: /*@C
218:    EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was 
219:    stopped.

221:    Not Collective

223:    Input Parameter:
224: .  eps - the eigensolver context

226:    Output Parameter:
227: .  reason - negative value indicates diverged, positive value converged
228:    (see EPSConvergedReason)

230:    Possible values for reason:
231: +  EPS_CONVERGED_TOL - converged up to tolerance
232: .  EPS_DIVERGED_ITS - required more than its to reach convergence
233: .  EPS_DIVERGED_BREAKDOWN - generic breakdown in method
234: -  EPS_DIVERGED_NONSYMMETRIC - The operator is nonsymmetric

236:    Level: intermediate

238:    Notes: Can only be called after the call to EPSSolve() is complete.

240: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
241: @*/
242: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
243: {
246:   *reason = eps->reason;
247:   return(0);
248: }

252: /*@
253:    EPSGetInvariantSubspace - Gets a basis of the computed invariant subspace.

255:    Not Collective

257:    Input Parameter:
258: .  eps - the eigensolver context
259:   
260:    Output Parameter:
261: .  v - an array of vectors

263:    Notes:
264:    This function should be called after EPSSolve() has finished.

266:    The user should provide in v an array of nconv vectors, where nconv is
267:    the value returned by EPSGetConverged().

269:    The vectors returned in v span an invariant subspace associated with the
270:    (nconv) computed eigenvalues. An invariant subspace X of A satisfies Ax 
271:    in X for all x in X (a similar definition applies for generalized 
272:    eigenproblems). 

274:    Level: intermediate

276: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve()
277: @*/
278: PetscErrorCode EPSGetInvariantSubspace(EPS eps, Vec *v)
279: {
280:   PetscErrorCode ierr;
281:   int            i;

286:   if (!eps->V) {
287:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
288:   }
289:   for (i=0;i<eps->nconv;i++) {
290:     VecCopy(eps->V[i],v[i]);
291:   }
292:   return(0);
293: }

297: /*@
298:    EPSGetEigenpair - Gets the i-th solution of the eigenproblem 
299:    as computed by EPSSolve(). The solution consists in both the eigenvalue
300:    and the eigenvector (if available).

302:    Not Collective

304:    Input Parameters:
305: +  eps - eigensolver context 
306: -  i   - index of the solution

308:    Output Parameters:
309: +  eigr - real part of eigenvalue
310: .  eigi - imaginary part of eigenvalue
311: .  Vr   - real part of eigenvector
312: -  Vi   - imaginary part of eigenvector

314:    Notes:
315:    If the eigenvalue is real, then eigi and Vi are set to zero. In the 
316:    complex case (e.g. with BOPT=O_complex) the eigenvalue is stored 
317:    directly in eigr (eigi is set to zero) and the eigenvector Vr (Vi is 
318:    set to zero).

320:    The index i should be a value between 0 and nconv (see EPSGetConverged()).
321:    Eigenpairs are indexed according to the ordering criterion established 
322:    with EPSSetWhichEigenpairs().

324:    Level: beginner

326: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), 
327:           EPSGetInvariantSubspace()
328: @*/
329: PetscErrorCode EPSGetEigenpair(EPS eps, int i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
330: {
331:   PetscErrorCode ierr;
332:   int            k;
333:   PetscScalar    zero = 0.0;
334: #ifndef PETSC_USE_COMPLEX
335:   PetscScalar    minus = -1.0;
336: #endif

340:   if (!eps->eigr || !eps->eigi || !eps->V) {
341:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
342:   }
343:   if (i<0 || i>=eps->nconv) {
344:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
345:   }
346:   if (!eps->evecsavailable && (Vr || Vi) ) {
347:     (*eps->ops->computevectors)(eps);
348:   }

350:   if (!eps->perm) k = i;
351:   else k = eps->perm[i];
352: #ifdef PETSC_USE_COMPLEX
353:   if (eigr) *eigr = eps->eigr[k];
354:   if (eigi) *eigi = 0;
355:   if (Vr) { VecCopy(eps->AV[k], Vr);  }
356:   if (Vi) { VecSet(&zero, Vi);  }
357: #else
358:   if (eigr) *eigr = eps->eigr[k];
359:   if (eigi) *eigi = eps->eigi[k];
360:   if (eps->eigi[k] > 0) { /* first value of conjugate pair */
361:     if (Vr) { VecCopy(eps->AV[k], Vr);  }
362:     if (Vi) { VecCopy(eps->AV[k+1], Vi);  }
363:   } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
364:     if (Vr) { VecCopy(eps->AV[k-1], Vr);  }
365:     if (Vi) {
366:       VecCopy(eps->AV[k], Vi);
367:       VecScale(&minus, Vi);
368:     }
369:   } else { /* real eigenvalue */
370:     if (Vr) { VecCopy(eps->AV[k], Vr);  }
371:     if (Vi) { VecSet(&zero, Vi);  }
372:   }
373: #endif
374: 
375:   return(0);
376: }

380: /*@
381:    EPSGetErrorEstimate - Returns the error bound associated to the i-th 
382:    approximate eigenpair.

384:    Not Collective

386:    Input Parameter:
387: +  eps - eigensolver context 
388: -  i   - index of eigenpair

390:    Output Parameter:
391: .  errest - the error estimate

393:    Level: advanced

395: .seealso: EPSComputeRelativeError()
396: @*/
397: PetscErrorCode EPSGetErrorEstimate(EPS eps, int i, PetscReal *errest)
398: {
401:   if (!eps->eigr || !eps->eigi) {
402:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
403:   }
404:   if (i<0 || i>=eps->nconv) {
405:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
406:   }
407:   if (eps->perm) i = eps->perm[i];
408:   if (errest) *errest = eps->errest[i];
409:   return(0);
410: }


415: /*@
416:    EPSComputeResidualNorm - Computes the residual norm associated with 
417:    the i-th converged approximate eigenpair.

419:    Collective on EPS

421:    Input Parameter:
422: .  eps - the eigensolver context
423: .  i   - the solution index

425:    Output Parameter:
426: .  norm - the residual norm, computed as ||Ax-kBx|| where k is the 
427:    eigenvalue and x is the eigenvector. 
428:    If k=0 then the residual norm is computed as ||Ax||.

430:    Notes:
431:    The index i should be a value between 0 and nconv (see EPSGetConverged()).
432:    Eigenpairs are indexed according to the ordering criterion established 
433:    with EPSSetWhichEigenpairs().

435:    Level: beginner

437: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
438: @*/
439: PetscErrorCode EPSComputeResidualNorm(EPS eps, int i, PetscReal *norm)
440: {
441:   PetscErrorCode ierr;
442:   Vec            u, v, w, xr, xi;
443:   Mat            A, B;
444:   PetscScalar    alpha, kr, ki;
445: #ifndef PETSC_USE_COMPLEX
446:   PetscReal      ni, nr;
447: #endif
448: 
451:   STGetOperators(eps->OP,&A,&B);
452:   VecDuplicate(eps->vec_initial,&u);
453:   VecDuplicate(eps->vec_initial,&v);
454:   VecDuplicate(eps->vec_initial,&w);
455:   VecDuplicate(eps->vec_initial,&xr);
456:   VecDuplicate(eps->vec_initial,&xi);
457:   EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);

459: #ifndef PETSC_USE_COMPLEX
460:   if (ki == 0 ||
461:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
462: #endif
463:     MatMult( A, xr, u );  /* u=A*x */
464:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
465:       if (eps->isgeneralized) { MatMult( B, xr, w );  }
466:       else { VecCopy( xr, w );  } /* w=B*x */
467:       alpha = -kr;
468:       VecAXPY( &alpha, w, u );  /* u=A*x-k*B*x */
469:     }
470:     VecNorm( u, NORM_2, norm);
471: #ifndef PETSC_USE_COMPLEX
472:   } else {
473:     MatMult( A, xr, u );  /* u=A*xr */
474:     if (eps->isgeneralized) { MatMult( B, xr, v );  }
475:     else { VecCopy( xr, v );  } /* v=B*xr */
476:     alpha = -kr;
477:     VecAXPY( &alpha, v, u );  /* u=A*xr-kr*B*xr */
478:     if (eps->isgeneralized) { MatMult( B, xi, w );  }
479:     else { VecCopy( xi, w );  } /* w=B*xi */
480:     alpha = ki;
481:     VecAXPY( &alpha, w, u );  /* u=A*xr-kr*B*xr+ki*B*xi */
482:     VecNorm( u, NORM_2, &nr );
483:     MatMult( A, xi, u );  /* u=A*xi */
484:     alpha = -kr;
485:     VecAXPY( &alpha, w, u );  /* u=A*xi-kr*B*xi */
486:     alpha = -ki;
487:     VecAXPY( &alpha, v, u );  /* u=A*xi-kr*B*xi-ki*B*xr */
488:     VecNorm( u, NORM_2, &ni );
489:     *norm = SlepcAbsEigenvalue( nr, ni );
490:   }
491: #endif

493:   VecDestroy(w);
494:   VecDestroy(v);
495:   VecDestroy(u);
496:   VecDestroy(xr);
497:   VecDestroy(xi);
498:   return(0);
499: }

503: /*@
504:    EPSComputeRelativeError - Computes the actual relative error associated 
505:    with the i-th converged approximate eigenpair.

507:    Collective on EPS

509:    Input Parameter:
510: .  eps - the eigensolver context
511: .  i   - the solution index

513:    Output Parameter:
514: .  error - the relative error, computed as ||Ax-kBx||/||kx|| where k is the 
515:    eigenvalue and x is the eigenvector. 
516:    If k=0 the relative error is computed as ||Ax||/||x||.

518:    Level: beginner

520: .seealso: EPSSolve(), EPSComputeResidualNorm()
521: @*/
522: PetscErrorCode EPSComputeRelativeError(EPS eps, int i, PetscReal *error)
523: {
524:   PetscErrorCode ierr;
525:   Vec            xr, xi;
526:   PetscScalar    kr, ki;
527:   PetscReal      norm, er;
528: #ifndef PETSC_USE_COMPLEX
529:   Vec            u;
530:   PetscScalar    alpha;
531:   PetscReal      ei;
532: #endif
533: 
536:   EPSComputeResidualNorm(eps,i,&norm);
537:   VecDuplicate(eps->vec_initial,&xr);
538:   VecDuplicate(eps->vec_initial,&xi);
539:   EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);

541: #ifndef PETSC_USE_COMPLEX
542:   if (ki == 0 ||
543:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
544: #endif
545:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
546:       VecScale(&kr, xr);
547:     }
548:     VecNorm(xr, NORM_2, &er);
549:     *error = norm / er;
550: #ifndef PETSC_USE_COMPLEX
551:   } else {
552:     VecDuplicate(xi, &u);
553:     VecCopy(xi, u);
554:     alpha = -ki;
555:     VecAXPBY(&kr, &alpha, xr, u);
556:     VecNorm(u, NORM_2, &er);
557:     VecAXPBY(&kr, &ki, xr, xi);
558:     VecNorm(xi, NORM_2, &ei);
559:     VecDestroy(u);
560:     *error = norm / SlepcAbsEigenvalue(er, ei);
561:   }
562: #endif    
563: 
564:   VecDestroy(xr);
565:   VecDestroy(xi);
566:   return(0);
567: }

571: /*@
572:    EPSReverseProjection - Compute the operation V=V*S, where the columns of
573:    V are m of the basis vectors of the EPS object and S is an mxm dense
574:    matrix.

576:    Collective on EPS

578:    Input Parameter:
579: +  eps - the eigenproblem solver context
580: .  V - basis vectors
581: .  S - pointer to the values of matrix S
582: .  k - starting column
583: .  m - dimension of matrix S
584: -  work - workarea of m vectors for intermediate results

586:    Level: developer

588: @*/
589: PetscErrorCode EPSReverseProjection(EPS eps,Vec* V,PetscScalar *S,int k,int m,Vec* work)
590: {
591:   PetscErrorCode ierr;
592:   int            i;
593:   PetscScalar    zero = 0.0;
594: 
596:   for (i=k;i<m;i++) {
597:     VecSet(&zero,work[i]);
598:     VecMAXPY(m,S+m*i,work[i],V);
599:   }
600:   for (i=k;i<m;i++) {
601:     VecCopy(work[i],V[i]);
602:   }
603:   return(0);
604: }

608: /*@
609:     EPSComputeExplicitOperator - Computes the explicit operator associated
610:     to the eigenvalue problem with the specified spectral transformation.  

612:     Collective on EPS

614:     Input Parameter:
615: .   eps - the eigenvalue solver context

617:     Output Parameter:
618: .   mat - the explicit operator

620:     Notes:
621:     This routine builds a matrix containing the explicit operator. For 
622:     example, in generalized problems with shift-and-invert spectral
623:     transformation the result would be matrix (A - s B)^-1 B.
624:     
625:     This computation is done by applying the operator to columns of the 
626:     identity matrix.

628:     Currently, this routine uses a dense matrix format when 1 processor
629:     is used and a sparse format otherwise.  This routine is costly in general,
630:     and is recommended for use only with relatively small systems.

632:     Level: advanced

634: .seealso: STApply()   
635: @*/
636: PetscErrorCode EPSComputeExplicitOperator(EPS eps,Mat *mat)
637: {
638:   PetscErrorCode ierr;
639:   Vec            in,out;
640:   int            i,M,m,size,*rows,start,end;
641:   MPI_Comm       comm;
642:   PetscScalar    *array,zero = 0.0,one = 1.0;

647:   comm = eps->comm;

649:   MPI_Comm_size(comm,&size);

651:   VecDuplicate(eps->vec_initial,&in);
652:   VecDuplicate(eps->vec_initial,&out);
653:   VecGetSize(in,&M);
654:   VecGetLocalSize(in,&m);
655:   VecGetOwnershipRange(in,&start,&end);
656:   PetscMalloc((m+1)*sizeof(int),&rows);
657:   for (i=0; i<m; i++) {rows[i] = start + i;}

659:   if (size == 1) {
660:     MatCreateSeqDense(comm,M,M,PETSC_NULL,mat);
661:   } else {
662:     MatCreateMPIAIJ(comm,m,m,M,M,0,0,0,0,mat);
663:   }
664: 
665:   for (i=0; i<M; i++) {

667:     VecSet(&zero,in);
668:     VecSetValues(in,1,&i,&one,INSERT_VALUES);
669:     VecAssemblyBegin(in);
670:     VecAssemblyEnd(in);

672:     STApply(eps->OP,in,out);
673: 
674:     VecGetArray(out,&array);
675:     MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
676:     VecRestoreArray(out,&array);

678:   }
679:   PetscFree(rows);
680:   VecDestroy(in);
681:   VecDestroy(out);
682:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
683:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
684:   return(0);
685: }