Actual source code: qepsolve.c

  1: /*
  2:       QEP routines related to the solution process.

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

  8:    This file is part of SLEPc.
  9:       
 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/qepimpl.h>       /*I "slepcqep.h" I*/

 28: /*@
 29:    QEPSolve - Solves the quadratic eigensystem.

 31:    Collective on QEP

 33:    Input Parameter:
 34: .  qep - eigensolver context obtained from QEPCreate()

 36:    Options Database:
 37: +   -qep_view - print information about the solver used
 38: .   -qep_view_binary - save the matrices to the default binary file
 39: -   -qep_plot_eigs - plot computed eigenvalues

 41:    Level: beginner

 43: .seealso: QEPCreate(), QEPSetUp(), QEPDestroy(), QEPSetTolerances() 
 44: @*/
 45: PetscErrorCode QEPSolve(QEP qep)
 46: {
 48:   PetscInt       i;
 49:   PetscReal      re,im;
 50:   PetscBool      flg,viewed=PETSC_FALSE;
 51:   PetscViewer    viewer;
 52:   PetscDraw      draw;
 53:   PetscDrawSP    drawsp;
 54:   char           filename[PETSC_MAX_PATH_LEN];
 55:   char           view[10];


 60:   flg = PETSC_FALSE;
 61:   PetscOptionsGetBool(((PetscObject)qep)->prefix,"-qep_view_binary",&flg,PETSC_NULL);
 62:   if (flg) {
 63:     MatView(qep->M,PETSC_VIEWER_BINARY_(((PetscObject)qep)->comm));
 64:     MatView(qep->C,PETSC_VIEWER_BINARY_(((PetscObject)qep)->comm));
 65:     MatView(qep->K,PETSC_VIEWER_BINARY_(((PetscObject)qep)->comm));
 66:   }

 68:   PetscOptionsGetString(((PetscObject)qep)->prefix,"-qep_view",view,10,&flg);
 69:   if (flg) {
 70:     PetscStrcmp(view,"before",&viewed);
 71:     if (viewed){
 72:       PetscViewer viewer;
 73:       PetscViewerASCIIGetStdout(((PetscObject)qep)->comm,&viewer);
 74:       QEPView(qep,viewer);
 75:     }
 76:   }

 78:   /* reset the convergence flag from the previous solves */
 79:   qep->reason = QEP_CONVERGED_ITERATING;

 81:   if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
 82:   if (!qep->setupcalled){ QEPSetUp(qep); }
 83:   qep->nconv = 0;
 84:   qep->its = 0;
 85:   for (i=0;i<qep->ncv;i++) qep->eigr[i]=qep->eigi[i]=qep->errest[i]=0.0;
 86:   QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,qep->ncv);

 88:   DSSetEigenvalueComparison(qep->ds,qep->which_func,qep->which_ctx);

 90:   PetscLogEventBegin(QEP_Solve,qep,0,0,0);
 91:   (*qep->ops->solve)(qep);
 92:   PetscLogEventEnd(QEP_Solve,qep,0,0,0);

 94:   if (!qep->reason) SETERRQ(((PetscObject)qep)->comm,1,"Internal error, solver returned without setting converged reason");

 96: #if !defined(PETSC_USE_COMPLEX)
 97:   /* reorder conjugate eigenvalues (positive imaginary first) */
 98:   for (i=0;i<qep->nconv-1;i++) {
 99:     if (qep->eigi[i] != 0) {
100:       if (qep->eigi[i] < 0) {
101:         qep->eigi[i] = -qep->eigi[i];
102:         qep->eigi[i+1] = -qep->eigi[i+1];
103:         VecScale(qep->V[i+1],-1.0);
104:       }
105:       i++;
106:     }
107:   }
108: #endif

110:   /* sort eigenvalues according to qep->which parameter */
111:   QEPSortEigenvalues(qep,qep->nconv,qep->eigr,qep->eigi,qep->perm);

113:   if (!viewed) {
114:     PetscOptionsGetString(((PetscObject)qep)->prefix,"-qep_view",filename,PETSC_MAX_PATH_LEN,&flg);
115:     if (flg && !PetscPreLoadingOn) {
116:       PetscViewerASCIIOpen(((PetscObject)qep)->comm,filename,&viewer);
117:       QEPView(qep,viewer);
118:       PetscViewerDestroy(&viewer);
119:     }
120:   }

122:   flg = PETSC_FALSE;
123:   PetscOptionsGetBool(((PetscObject)qep)->prefix,"-qep_plot_eigs",&flg,PETSC_NULL);
124:   if (flg) {
125:     PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
126:                              PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
127:     PetscViewerDrawGetDraw(viewer,0,&draw);
128:     PetscDrawSPCreate(draw,1,&drawsp);
129:     for(i=0;i<qep->nconv;i++) {
130: #if defined(PETSC_USE_COMPLEX)
131:       re = PetscRealPart(qep->eigr[i]);
132:       im = PetscImaginaryPart(qep->eigi[i]);
133: #else
134:       re = qep->eigr[i];
135:       im = qep->eigi[i];
136: #endif
137:       PetscDrawSPAddPoint(drawsp,&re,&im);
138:     }
139:     PetscDrawSPDraw(drawsp);
140:     PetscDrawSPDestroy(&drawsp);
141:     PetscViewerDestroy(&viewer);
142:   }

144:   /* Remove the initial subspace */
145:   qep->nini = 0;
146:   return(0);
147: }

151: /*@
152:    QEPGetIterationNumber - Gets the current iteration number. If the 
153:    call to QEPSolve() is complete, then it returns the number of iterations 
154:    carried out by the solution method.
155:  
156:    Not Collective

158:    Input Parameter:
159: .  qep - the quadratic eigensolver context

161:    Output Parameter:
162: .  its - number of iterations

164:    Level: intermediate

166:    Note:
167:    During the i-th iteration this call returns i-1. If QEPSolve() is 
168:    complete, then parameter "its" contains either the iteration number at
169:    which convergence was successfully reached, or failure was detected.  
170:    Call QEPGetConvergedReason() to determine if the solver converged or 
171:    failed and why.

173: .seealso: QEPGetConvergedReason(), QEPSetTolerances()
174: @*/
175: PetscErrorCode QEPGetIterationNumber(QEP qep,PetscInt *its)
176: {
180:   *its = qep->its;
181:   return(0);
182: }

186: /*@
187:    QEPGetConverged - Gets the number of converged eigenpairs.

189:    Not Collective

191:    Input Parameter:
192: .  qep - the quadratic eigensolver context
193:   
194:    Output Parameter:
195: .  nconv - number of converged eigenpairs 

197:    Note:
198:    This function should be called after QEPSolve() has finished.

200:    Level: beginner

202: .seealso: QEPSetDimensions(), QEPSolve()
203: @*/
204: PetscErrorCode QEPGetConverged(QEP qep,PetscInt *nconv)
205: {
209:   *nconv = qep->nconv;
210:   return(0);
211: }

215: /*@C
216:    QEPGetConvergedReason - Gets the reason why the QEPSolve() iteration was 
217:    stopped.

219:    Not Collective

221:    Input Parameter:
222: .  qep - the quadratic eigensolver context

224:    Output Parameter:
225: .  reason - negative value indicates diverged, positive value converged

227:    Possible values for reason:
228: +  QEP_CONVERGED_TOL - converged up to tolerance
229: .  QEP_DIVERGED_ITS - required more than its to reach convergence
230: -  QEP_DIVERGED_BREAKDOWN - generic breakdown in method

232:    Note:
233:    Can only be called after the call to QEPSolve() is complete.

235:    Level: intermediate

237: .seealso: QEPSetTolerances(), QEPSolve(), QEPConvergedReason
238: @*/
239: PetscErrorCode QEPGetConvergedReason(QEP qep,QEPConvergedReason *reason)
240: {
244:   *reason = qep->reason;
245:   return(0);
246: }

250: /*@
251:    QEPGetEigenpair - Gets the i-th solution of the eigenproblem as computed by 
252:    QEPSolve(). The solution consists in both the eigenvalue and the eigenvector.

254:    Logically Collective on EPS

256:    Input Parameters:
257: +  qep - quadratic eigensolver context 
258: -  i   - index of the solution

260:    Output Parameters:
261: +  eigr - real part of eigenvalue
262: .  eigi - imaginary part of eigenvalue
263: .  Vr   - real part of eigenvector
264: -  Vi   - imaginary part of eigenvector

266:    Notes:
267:    If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is 
268:    configured with complex scalars the eigenvalue is stored 
269:    directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is 
270:    set to zero).

272:    The index i should be a value between 0 and nconv-1 (see QEPGetConverged()).
273:    Eigenpairs are indexed according to the ordering criterion established 
274:    with QEPSetWhichEigenpairs().

276:    Level: beginner

278: .seealso: QEPSolve(), QEPGetConverged(), QEPSetWhichEigenpairs()
279: @*/
280: PetscErrorCode QEPGetEigenpair(QEP qep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
281: {
282:   PetscInt       k;

290:   if (!qep->eigr || !qep->eigi || !qep->V) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONGSTATE,"QEPSolve must be called first");
291:   if (i<0 || i>=qep->nconv) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");

293:   if (!qep->perm) k = i;
294:   else k = qep->perm[i];

296:   /* eigenvalue */
297: #if defined(PETSC_USE_COMPLEX)
298:   if (eigr) *eigr = qep->eigr[k];
299:   if (eigi) *eigi = 0;
300: #else
301:   if (eigr) *eigr = qep->eigr[k];
302:   if (eigi) *eigi = qep->eigi[k];
303: #endif
304: 
305:   /* eigenvector */
306: #if defined(PETSC_USE_COMPLEX)
307:   if (Vr) { VecCopy(qep->V[k],Vr); }
308:   if (Vi) { VecSet(Vi,0.0); }
309: #else
310:   if (qep->eigi[k]>0) { /* first value of conjugate pair */
311:     if (Vr) { VecCopy(qep->V[k],Vr); }
312:     if (Vi) { VecCopy(qep->V[k+1],Vi); }
313:   } else if (qep->eigi[k]<0) { /* second value of conjugate pair */
314:     if (Vr) { VecCopy(qep->V[k-1],Vr); }
315:     if (Vi) {
316:       VecCopy(qep->V[k],Vi);
317:       VecScale(Vi,-1.0);
318:     }
319:   } else { /* real eigenvalue */
320:     if (Vr) { VecCopy(qep->V[k],Vr); }
321:     if (Vi) { VecSet(Vi,0.0); }
322:   }
323: #endif
324:   return(0);
325: }

329: /*@
330:    QEPGetErrorEstimate - Returns the error estimate associated to the i-th 
331:    computed eigenpair.

333:    Not Collective

335:    Input Parameter:
336: +  qep - quadratic eigensolver context 
337: -  i   - index of eigenpair

339:    Output Parameter:
340: .  errest - the error estimate

342:    Notes:
343:    This is the error estimate used internally by the eigensolver. The actual
344:    error bound can be computed with QEPComputeRelativeError(). See also the users
345:    manual for details.

347:    Level: advanced

349: .seealso: QEPComputeRelativeError()
350: @*/
351: PetscErrorCode QEPGetErrorEstimate(QEP qep,PetscInt i,PetscReal *errest)
352: {
356:   if (!qep->eigr || !qep->eigi) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"QEPSolve must be called first");
357:   if (i<0 || i>=qep->nconv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
358:   if (qep->perm) i = qep->perm[i];
359:   if (errest) *errest = qep->errest[i];
360:   return(0);
361: }

365: /*
366:    QEPComputeResidualNorm_Private - Computes the norm of the residual vector 
367:    associated with an eigenpair.
368: */
369: PetscErrorCode QEPComputeResidualNorm_Private(QEP qep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *norm)
370: {
372:   Vec            u,w;
373:   Mat            M=qep->M,C=qep->C,K=qep->K;
374: #if !defined(PETSC_USE_COMPLEX)
375:   Vec            v,y,z;
376:   PetscReal      ni,nr;
377:   PetscScalar    a1,a2;
378: #endif
379: 
381:   VecDuplicate(qep->V[0],&u);
382:   VecDuplicate(u,&w);
383: 
384: #if !defined(PETSC_USE_COMPLEX)
385:   if (ki == 0 ||
386:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
387: #endif
388:     MatMult(K,xr,u);                 /* u=K*x */
389:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
390:       MatMult(C,xr,w);               /* w=C*x */
391:       VecAXPY(u,kr,w);               /* u=l*C*x+K*x */
392:       MatMult(M,xr,w);               /* w=M*x */
393:       VecAXPY(u,kr*kr,w);            /* u=l^2*M*x+l*C*x+K*x */
394:     }
395:     VecNorm(u,NORM_2,norm);
396: #if !defined(PETSC_USE_COMPLEX)
397:   } else {
398:     VecDuplicate(u,&v);
399:     VecDuplicate(u,&y);
400:     VecDuplicate(u,&z);
401:     a1 = kr*kr-ki*ki;
402:     a2 = 2.0*kr*ki;
403:     MatMult(K,xr,u);           /* u=K*xr */
404:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
405:       MatMult(C,xr,v);         /* v=C*xr */
406:       MatMult(C,xi,w);         /* w=C*xi */
407:       MatMult(M,xr,y);         /* y=M*xr */
408:       MatMult(M,xi,z);         /* z=M*xi */
409:       VecAXPY(u,kr,v);         /* u=kr*C*xr+K*xr */
410:       VecAXPY(u,-ki,w);        /* u=kr*C*xr-ki*C*xi+K*xr */
411:       VecAXPY(u,a1,y);         /* u=a1*M*xr+kr*C*xr-ki*C*xi+K*xr */
412:       VecAXPY(u,-a2,z);        /* u=a1*M*xr-a2*M*xi+kr*C*xr-ki*C*xi+K*xr */
413:     }
414:     VecNorm(u,NORM_2,&nr);
415:     MatMult(K,xi,u);         /* u=K*xi */
416:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
417:       VecAXPY(u,kr,w);         /* u=kr*C*xi+K*xi */
418:       VecAXPY(u,ki,v);         /* u=kr*C*xi+ki*C*xi+K*xi */
419:       VecAXPY(u,a1,z);         /* u=a1*M*xi+kr*C*xi+ki*C*xi+K*xi */
420:       VecAXPY(u,a2,y);         /* u=a1*M*xi+a2*M*ki+kr*C*xi+ki*C*xi+K*xi */
421:     }
422:     VecNorm(u,NORM_2,&ni);
423:     *norm = SlepcAbsEigenvalue(nr,ni);
424:     VecDestroy(&v);
425:     VecDestroy(&y);
426:     VecDestroy(&z);
427:   }
428: #endif

430:   VecDestroy(&w);
431:   VecDestroy(&u);
432:   return(0);
433: }

437: /*@
438:    QEPComputeResidualNorm - Computes the norm of the residual vector associated with 
439:    the i-th computed eigenpair.

441:    Collective on QEP

443:    Input Parameter:
444: .  qep - the quadratic eigensolver context
445: .  i   - the solution index

447:    Output Parameter:
448: .  norm - the residual norm, computed as ||(l^2*M+l*C+K)x||_2 where l is the 
449:    eigenvalue and x is the eigenvector. 
450:    If l=0 then the residual norm is computed as ||Kx||_2.

452:    Notes:
453:    The index i should be a value between 0 and nconv-1 (see QEPGetConverged()).
454:    Eigenpairs are indexed according to the ordering criterion established 
455:    with QEPSetWhichEigenpairs().

457:    Level: beginner

459: .seealso: QEPSolve(), QEPGetConverged(), QEPSetWhichEigenpairs()
460: @*/
461: PetscErrorCode QEPComputeResidualNorm(QEP qep,PetscInt i,PetscReal *norm)
462: {
464:   Vec            xr,xi;
465:   PetscScalar    kr,ki;
466: 
471:   VecDuplicate(qep->V[0],&xr);
472:   VecDuplicate(qep->V[0],&xi);
473:   QEPGetEigenpair(qep,i,&kr,&ki,xr,xi);
474:   QEPComputeResidualNorm_Private(qep,kr,ki,xr,xi,norm);
475:   VecDestroy(&xr);
476:   VecDestroy(&xi);
477:   return(0);
478: }

482: /*
483:    QEPComputeRelativeError_Private - Computes the relative error bound 
484:    associated with an eigenpair.
485: */
486: PetscErrorCode QEPComputeRelativeError_Private(QEP qep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *error)
487: {
489:   PetscReal      norm,er;
490: #if !defined(PETSC_USE_COMPLEX)
491:   PetscReal      ei;
492: #endif
493: 
495:   QEPComputeResidualNorm_Private(qep,kr,ki,xr,xi,&norm);
496: #if !defined(PETSC_USE_COMPLEX)
497:   if (ki == 0 ||
498:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
499: #endif
500:     VecNorm(xr,NORM_2,&er);
501:     if (PetscAbsScalar(kr) > norm) {
502:       *error = norm/(PetscAbsScalar(kr)*er);
503:     } else {
504:       *error = norm/er;
505:     }
506: #if !defined(PETSC_USE_COMPLEX)
507:   } else {
508:     VecNorm(xr,NORM_2,&er);
509:     VecNorm(xi,NORM_2,&ei);
510:     if (SlepcAbsEigenvalue(kr,ki) > norm) {
511:       *error = norm/(SlepcAbsEigenvalue(kr,ki)*SlepcAbsEigenvalue(er,ei));
512:     } else {
513:       *error = norm/SlepcAbsEigenvalue(er,ei);
514:     }
515:   }
516: #endif    
517:   return(0);
518: }

522: /*@
523:    QEPComputeRelativeError - Computes the relative error bound associated 
524:    with the i-th computed eigenpair.

526:    Collective on QEP

528:    Input Parameter:
529: .  qep - the quadratic eigensolver context
530: .  i   - the solution index

532:    Output Parameter:
533: .  error - the relative error bound, computed as ||(l^2*M+l*C+K)x||_2/||lx||_2 where 
534:    l is the eigenvalue and x is the eigenvector. 
535:    If l=0 the relative error is computed as ||Kx||_2/||x||_2.

537:    Level: beginner

539: .seealso: QEPSolve(), QEPComputeResidualNorm(), QEPGetErrorEstimate()
540: @*/
541: PetscErrorCode QEPComputeRelativeError(QEP qep,PetscInt i,PetscReal *error)
542: {
544:   Vec            xr,xi;
545:   PetscScalar    kr,ki;
546: 
551:   VecDuplicate(qep->V[0],&xr);
552:   VecDuplicate(qep->V[0],&xi);
553:   QEPGetEigenpair(qep,i,&kr,&ki,xr,xi);
554:   QEPComputeRelativeError_Private(qep,kr,ki,xr,xi,error);
555:   VecDestroy(&xr);
556:   VecDestroy(&xi);
557:   return(0);
558: }

562: /*@
563:    QEPSortEigenvalues - Sorts a list of eigenvalues according to the criterion 
564:    specified via QEPSetWhichEigenpairs().

566:    Not Collective

568:    Input Parameters:
569: +  qep   - the quadratic eigensolver context
570: .  n     - number of eigenvalues in the list
571: .  eigr  - pointer to the array containing the eigenvalues
572: -  eigi  - imaginary part of the eigenvalues (only when using real numbers)

574:    Output Parameter:
575: .  perm  - resulting permutation

577:    Note:
578:    The result is a list of indices in the original eigenvalue array 
579:    corresponding to the first nev eigenvalues sorted in the specified
580:    criterion.

582:    Level: developer

584: .seealso: QEPSetWhichEigenpairs()
585: @*/
586: PetscErrorCode QEPSortEigenvalues(QEP qep,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm)
587: {
589:   PetscScalar    re,im;
590:   PetscInt       i,j,result,tmp;

597:   for (i=0; i<n; i++) { perm[i] = i; }
598:   /* insertion sort */
599:   for (i=n-1; i>=0; i--) {
600:     re = eigr[perm[i]];
601:     im = eigi[perm[i]];
602:     j = i + 1;
603: #if !defined(PETSC_USE_COMPLEX)
604:     if (im != 0) {
605:       /* complex eigenvalue */
606:       i--;
607:       im = eigi[perm[i]];
608:     }
609: #endif
610:     while (j<n) {
611:       QEPCompareEigenvalues(qep,re,im,eigr[perm[j]],eigi[perm[j]],&result);
612:       if (result < 0) break;
613: #if !defined(PETSC_USE_COMPLEX)
614:       /* keep together every complex conjugated eigenpair */
615:       if (im == 0) {
616:         if (eigi[perm[j]] == 0) {
617: #endif
618:           tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = tmp;
619:           j++;
620: #if !defined(PETSC_USE_COMPLEX)
621:         } else {
622:           tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp;
623:           j+=2;
624:         }
625:       } else {
626:         if (eigi[perm[j]] == 0) {
627:           tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = perm[j-1]; perm[j-1] = tmp;
628:           j++;
629:         } else {
630:           tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = tmp;
631:           tmp = perm[j-1]; perm[j-1] = perm[j+1]; perm[j+1] = tmp;
632:           j+=2;
633:         }
634:       }
635: #endif
636:     }
637:   }
638:   return(0);
639: }

643: /*@
644:    QEPCompareEigenvalues - Compares two (possibly complex) eigenvalues according
645:    to a certain criterion.

647:    Not Collective

649:    Input Parameters:
650: +  qep    - the quadratic eigensolver context
651: .  ar     - real part of the 1st eigenvalue
652: .  ai     - imaginary part of the 1st eigenvalue
653: .  br     - real part of the 2nd eigenvalue
654: -  bi     - imaginary part of the 2nd eigenvalue

656:    Output Parameter:
657: .  res    - result of comparison

659:    Notes:
660:    Returns an integer less than, equal to, or greater than zero if the first
661:    eigenvalue is considered to be respectively less than, equal to, or greater
662:    than the second one.

664:    The criterion of comparison is related to the 'which' parameter set with
665:    QEPSetWhichEigenpairs().

667:    Level: developer

669: .seealso: QEPSortEigenvalues(), QEPSetWhichEigenpairs()
670: @*/
671: PetscErrorCode QEPCompareEigenvalues(QEP qep,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result)
672: {

678:   if (!qep->which_func) SETERRQ(PETSC_COMM_SELF,1,"Undefined eigenvalue comparison function");
679:   (*qep->which_func)(ar,ai,br,bi,result,qep->which_ctx);
680:   return(0);
681: }

685: /*@
686:    QEPGetOperationCounters - Gets the total number of matrix-vector products, dot 
687:    products, and linear solve iterations used by the QEP object during the last
688:    QEPSolve() call.

690:    Not Collective

692:    Input Parameter:
693: .  qep - quadratic eigensolver context

695:    Output Parameter:
696: +  matvecs - number of matrix-vector product operations
697: .  dots    - number of dot product operations
698: -  lits    - number of linear iterations

700:    Notes:
701:    These counters are reset to zero at each successive call to QEPSolve().

703:    Level: intermediate

705: @*/
706: PetscErrorCode QEPGetOperationCounters(QEP qep,PetscInt* matvecs,PetscInt* dots,PetscInt* lits)
707: {
709: 
712:   if (matvecs) *matvecs = qep->matvecs;
713:   if (dots) {
714:     if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
715:     IPGetOperationCounters(qep->ip,dots);
716:   }
717:   if (lits) *lits = qep->linits;
718:   return(0);
719: }