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-2011, 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 <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;
51: PetscViewer viewer;
52: PetscDraw draw;
53: PetscDrawSP drawsp;
54: char filename[PETSC_MAX_PATH_LEN];
59: flg = PETSC_FALSE;
60: PetscOptionsGetBool(((PetscObject)qep)->prefix,"-qep_view_binary",&flg,PETSC_NULL);
61: if (flg) {
62: MatView(qep->M,PETSC_VIEWER_BINARY_(((PetscObject)qep)->comm));
63: MatView(qep->C,PETSC_VIEWER_BINARY_(((PetscObject)qep)->comm));
64: MatView(qep->K,PETSC_VIEWER_BINARY_(((PetscObject)qep)->comm));
65: }
67: /* reset the convergence flag from the previous solves */
68: qep->reason = QEP_CONVERGED_ITERATING;
70: if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
71: if (!qep->setupcalled){ QEPSetUp(qep); }
72: qep->nconv = 0;
73: qep->its = 0;
74: for (i=0;i<qep->ncv;i++) qep->eigr[i]=qep->eigi[i]=qep->errest[i]=0.0;
75: QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,qep->ncv);
77: PetscLogEventBegin(QEP_Solve,qep,0,0,0);
78: (*qep->ops->solve)(qep);
79: PetscLogEventEnd(QEP_Solve,qep,0,0,0);
81: if (!qep->reason) {
82: SETERRQ(((PetscObject)qep)->comm,1,"Internal error, solver returned without setting converged reason");
83: }
85: #if !defined(PETSC_USE_COMPLEX)
86: /* reorder conjugate eigenvalues (positive imaginary first) */
87: for (i=0;i<qep->nconv-1;i++) {
88: if (qep->eigi[i] != 0) {
89: if (qep->eigi[i] < 0) {
90: qep->eigi[i] = -qep->eigi[i];
91: qep->eigi[i+1] = -qep->eigi[i+1];
92: VecScale(qep->V[i+1],-1.0);
93: }
94: i++;
95: }
96: }
97: #endif
99: /* sort eigenvalues according to qep->which parameter */
100: PetscFree(qep->perm);
101: if (qep->nconv > 0) {
102: PetscMalloc(sizeof(PetscInt)*qep->nconv,&qep->perm);
103: QEPSortEigenvalues(qep,qep->nconv,qep->eigr,qep->eigi,qep->perm);
104: }
106: PetscOptionsGetString(((PetscObject)qep)->prefix,"-qep_view",filename,PETSC_MAX_PATH_LEN,&flg);
107: if (flg && !PetscPreLoadingOn) {
108: PetscViewerASCIIOpen(((PetscObject)qep)->comm,filename,&viewer);
109: QEPView(qep,viewer);
110: PetscViewerDestroy(&viewer);
111: }
113: flg = PETSC_FALSE;
114: PetscOptionsGetBool(((PetscObject)qep)->prefix,"-qep_plot_eigs",&flg,PETSC_NULL);
115: if (flg) {
116: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
117: PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
118: PetscViewerDrawGetDraw(viewer,0,&draw);
119: PetscDrawSPCreate(draw,1,&drawsp);
120: for(i=0;i<qep->nconv;i++) {
121: #if defined(PETSC_USE_COMPLEX)
122: re = PetscRealPart(qep->eigr[i]);
123: im = PetscImaginaryPart(qep->eigi[i]);
124: #else
125: re = qep->eigr[i];
126: im = qep->eigi[i];
127: #endif
128: PetscDrawSPAddPoint(drawsp,&re,&im);
129: }
130: PetscDrawSPDraw(drawsp);
131: PetscDrawSPDestroy(&drawsp);
132: PetscViewerDestroy(&viewer);
133: }
135: /* Remove the initial subspace */
136: qep->nini = 0;
137: return(0);
138: }
142: /*@
143: QEPGetIterationNumber - Gets the current iteration number. If the
144: call to QEPSolve() is complete, then it returns the number of iterations
145: carried out by the solution method.
146:
147: Not Collective
149: Input Parameter:
150: . qep - the quadratic eigensolver context
152: Output Parameter:
153: . its - number of iterations
155: Level: intermediate
157: Note:
158: During the i-th iteration this call returns i-1. If QEPSolve() is
159: complete, then parameter "its" contains either the iteration number at
160: which convergence was successfully reached, or failure was detected.
161: Call QEPGetConvergedReason() to determine if the solver converged or
162: failed and why.
164: .seealso: QEPGetConvergedReason(), QEPSetTolerances()
165: @*/
166: PetscErrorCode QEPGetIterationNumber(QEP qep,PetscInt *its)
167: {
171: *its = qep->its;
172: return(0);
173: }
177: /*@
178: QEPGetConverged - Gets the number of converged eigenpairs.
180: Not Collective
182: Input Parameter:
183: . qep - the quadratic eigensolver context
184:
185: Output Parameter:
186: . nconv - number of converged eigenpairs
188: Note:
189: This function should be called after QEPSolve() has finished.
191: Level: beginner
193: .seealso: QEPSetDimensions(), QEPSolve()
194: @*/
195: PetscErrorCode QEPGetConverged(QEP qep,PetscInt *nconv)
196: {
200: *nconv = qep->nconv;
201: return(0);
202: }
206: /*@C
207: QEPGetConvergedReason - Gets the reason why the QEPSolve() iteration was
208: stopped.
210: Not Collective
212: Input Parameter:
213: . qep - the quadratic eigensolver context
215: Output Parameter:
216: . reason - negative value indicates diverged, positive value converged
218: Possible values for reason:
219: + QEP_CONVERGED_TOL - converged up to tolerance
220: . QEP_DIVERGED_ITS - required more than its to reach convergence
221: - QEP_DIVERGED_BREAKDOWN - generic breakdown in method
223: Note:
224: Can only be called after the call to QEPSolve() is complete.
226: Level: intermediate
228: .seealso: QEPSetTolerances(), QEPSolve(), QEPConvergedReason
229: @*/
230: PetscErrorCode QEPGetConvergedReason(QEP qep,QEPConvergedReason *reason)
231: {
235: *reason = qep->reason;
236: return(0);
237: }
241: /*@
242: QEPGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
243: QEPSolve(). The solution consists in both the eigenvalue and the eigenvector.
245: Not Collective, but vectors are shared by all processors that share the QEP
247: Input Parameters:
248: + qep - quadratic eigensolver context
249: - i - index of the solution
251: Output Parameters:
252: + eigr - real part of eigenvalue
253: . eigi - imaginary part of eigenvalue
254: . Vr - real part of eigenvector
255: - Vi - imaginary part of eigenvector
257: Notes:
258: If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
259: configured with complex scalars the eigenvalue is stored
260: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
261: set to zero).
263: The index i should be a value between 0 and nconv-1 (see QEPGetConverged()).
264: Eigenpairs are indexed according to the ordering criterion established
265: with QEPSetWhichEigenpairs().
267: Level: beginner
269: .seealso: QEPSolve(), QEPGetConverged(), QEPSetWhichEigenpairs()
270: @*/
271: PetscErrorCode QEPGetEigenpair(QEP qep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
272: {
273: PetscInt k;
280: if (!qep->eigr || !qep->eigi || !qep->V) {
281: SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONGSTATE,"QEPSolve must be called first");
282: }
283: if (i<0 || i>=qep->nconv) {
284: SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
285: }
287: if (!qep->perm) k = i;
288: else k = qep->perm[i];
290: /* eigenvalue */
291: #if defined(PETSC_USE_COMPLEX)
292: if (eigr) *eigr = qep->eigr[k];
293: if (eigi) *eigi = 0;
294: #else
295: if (eigr) *eigr = qep->eigr[k];
296: if (eigi) *eigi = qep->eigi[k];
297: #endif
298:
299: /* eigenvector */
300: #if defined(PETSC_USE_COMPLEX)
301: if (Vr) { VecCopy(qep->V[k],Vr); }
302: if (Vi) { VecSet(Vi,0.0); }
303: #else
304: if (qep->eigi[k]>0) { /* first value of conjugate pair */
305: if (Vr) { VecCopy(qep->V[k],Vr); }
306: if (Vi) { VecCopy(qep->V[k+1],Vi); }
307: } else if (qep->eigi[k]<0) { /* second value of conjugate pair */
308: if (Vr) { VecCopy(qep->V[k-1],Vr); }
309: if (Vi) {
310: VecCopy(qep->V[k],Vi);
311: VecScale(Vi,-1.0);
312: }
313: } else { /* real eigenvalue */
314: if (Vr) { VecCopy(qep->V[k],Vr); }
315: if (Vi) { VecSet(Vi,0.0); }
316: }
317: #endif
318: return(0);
319: }
323: /*@
324: QEPGetErrorEstimate - Returns the error estimate associated to the i-th
325: computed eigenpair.
327: Not Collective
329: Input Parameter:
330: + qep - quadratic eigensolver context
331: - i - index of eigenpair
333: Output Parameter:
334: . errest - the error estimate
336: Notes:
337: This is the error estimate used internally by the eigensolver. The actual
338: error bound can be computed with QEPComputeRelativeError(). See also the users
339: manual for details.
341: Level: advanced
343: .seealso: QEPComputeRelativeError()
344: @*/
345: PetscErrorCode QEPGetErrorEstimate(QEP qep,PetscInt i,PetscReal *errest)
346: {
350: if (!qep->eigr || !qep->eigi) {
351: SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONGSTATE,"QEPSolve must be called first");
352: }
353: if (i<0 || i>=qep->nconv) {
354: SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
355: }
356: if (qep->perm) i = qep->perm[i];
357: if (errest) *errest = qep->errest[i];
358: return(0);
359: }
363: /*
364: QEPComputeResidualNorm_Private - Computes the norm of the residual vector
365: associated with an eigenpair.
366: */
367: PetscErrorCode QEPComputeResidualNorm_Private(QEP qep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *norm)
368: {
370: Vec u,w;
371: Mat M=qep->M,C=qep->C,K=qep->K;
372: #if !defined(PETSC_USE_COMPLEX)
373: Vec v,y,z;
374: PetscReal ni,nr;
375: PetscScalar a1,a2;
376: #endif
377:
379: VecDuplicate(qep->V[0],&u);
380: VecDuplicate(u,&w);
381:
382: #if !defined(PETSC_USE_COMPLEX)
383: if (ki == 0 ||
384: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
385: #endif
386: MatMult(K,xr,u); /* u=K*x */
387: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
388: MatMult(C,xr,w); /* w=C*x */
389: VecAXPY(u,kr,w); /* u=l*C*x+K*x */
390: MatMult(M,xr,w); /* w=M*x */
391: VecAXPY(u,kr*kr,w); /* u=l^2*M*x+l*C*x+K*x */
392: }
393: VecNorm(u,NORM_2,norm);
394: #if !defined(PETSC_USE_COMPLEX)
395: } else {
396: VecDuplicate(u,&v);
397: VecDuplicate(u,&y);
398: VecDuplicate(u,&z);
399: a1 = kr*kr-ki*ki;
400: a2 = 2.0*kr*ki;
401: MatMult(K,xr,u); /* u=K*xr */
402: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
403: MatMult(C,xr,v); /* v=C*xr */
404: MatMult(C,xi,w); /* w=C*xi */
405: MatMult(M,xr,y); /* y=M*xr */
406: MatMult(M,xi,z); /* z=M*xi */
407: VecAXPY(u,kr,v); /* u=kr*C*xr+K*xr */
408: VecAXPY(u,-ki,w); /* u=kr*C*xr-ki*C*xi+K*xr */
409: VecAXPY(u,a1,y); /* u=a1*M*xr+kr*C*xr-ki*C*xi+K*xr */
410: VecAXPY(u,-a2,z); /* u=a1*M*xr-a2*M*xi+kr*C*xr-ki*C*xi+K*xr */
411: }
412: VecNorm(u,NORM_2,&nr);
413: MatMult(K,xi,u); /* u=K*xi */
414: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
415: VecAXPY(u,kr,w); /* u=kr*C*xi+K*xi */
416: VecAXPY(u,ki,v); /* u=kr*C*xi+ki*C*xi+K*xi */
417: VecAXPY(u,a1,z); /* u=a1*M*xi+kr*C*xi+ki*C*xi+K*xi */
418: VecAXPY(u,a2,y); /* u=a1*M*xi+a2*M*ki+kr*C*xi+ki*C*xi+K*xi */
419: }
420: VecNorm(u,NORM_2,&ni);
421: *norm = SlepcAbsEigenvalue(nr,ni);
422: VecDestroy(&v);
423: VecDestroy(&y);
424: VecDestroy(&z);
425: }
426: #endif
428: VecDestroy(&w);
429: VecDestroy(&u);
430: return(0);
431: }
435: /*@
436: QEPComputeResidualNorm - Computes the norm of the residual vector associated with
437: the i-th computed eigenpair.
439: Collective on QEP
441: Input Parameter:
442: . qep - the quadratic eigensolver context
443: . i - the solution index
445: Output Parameter:
446: . norm - the residual norm, computed as ||(l^2*M+l*C+K)x||_2 where l is the
447: eigenvalue and x is the eigenvector.
448: If l=0 then the residual norm is computed as ||Kx||_2.
450: Notes:
451: The index i should be a value between 0 and nconv-1 (see QEPGetConverged()).
452: Eigenpairs are indexed according to the ordering criterion established
453: with QEPSetWhichEigenpairs().
455: Level: beginner
457: .seealso: QEPSolve(), QEPGetConverged(), QEPSetWhichEigenpairs()
458: @*/
459: PetscErrorCode QEPComputeResidualNorm(QEP qep,PetscInt i,PetscReal *norm)
460: {
462: Vec xr,xi;
463: PetscScalar kr,ki;
464:
469: VecDuplicate(qep->V[0],&xr);
470: VecDuplicate(qep->V[0],&xi);
471: QEPGetEigenpair(qep,i,&kr,&ki,xr,xi);
472: QEPComputeResidualNorm_Private(qep,kr,ki,xr,xi,norm);
473: VecDestroy(&xr);
474: VecDestroy(&xi);
475: return(0);
476: }
480: /*
481: QEPComputeRelativeError_Private - Computes the relative error bound
482: associated with an eigenpair.
483: */
484: PetscErrorCode QEPComputeRelativeError_Private(QEP qep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *error)
485: {
487: PetscReal norm,er;
488: #if !defined(PETSC_USE_COMPLEX)
489: PetscReal ei;
490: #endif
491:
493: QEPComputeResidualNorm_Private(qep,kr,ki,xr,xi,&norm);
494: #if !defined(PETSC_USE_COMPLEX)
495: if (ki == 0 ||
496: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
497: #endif
498: VecNorm(xr,NORM_2,&er);
499: if (PetscAbsScalar(kr) > norm) {
500: *error = norm/(PetscAbsScalar(kr)*er);
501: } else {
502: *error = norm/er;
503: }
504: #if !defined(PETSC_USE_COMPLEX)
505: } else {
506: VecNorm(xr,NORM_2,&er);
507: VecNorm(xi,NORM_2,&ei);
508: if (SlepcAbsEigenvalue(kr,ki) > norm) {
509: *error = norm/(SlepcAbsEigenvalue(kr,ki)*SlepcAbsEigenvalue(er,ei));
510: } else {
511: *error = norm/SlepcAbsEigenvalue(er,ei);
512: }
513: }
514: #endif
515: return(0);
516: }
520: /*@
521: QEPComputeRelativeError - Computes the relative error bound associated
522: with the i-th computed eigenpair.
524: Collective on QEP
526: Input Parameter:
527: . qep - the quadratic eigensolver context
528: . i - the solution index
530: Output Parameter:
531: . error - the relative error bound, computed as ||(l^2*M+l*C+K)x||_2/||lx||_2 where
532: l is the eigenvalue and x is the eigenvector.
533: If l=0 the relative error is computed as ||Kx||_2/||x||_2.
535: Level: beginner
537: .seealso: QEPSolve(), QEPComputeResidualNorm(), QEPGetErrorEstimate()
538: @*/
539: PetscErrorCode QEPComputeRelativeError(QEP qep,PetscInt i,PetscReal *error)
540: {
542: Vec xr,xi;
543: PetscScalar kr,ki;
544:
549: VecDuplicate(qep->V[0],&xr);
550: VecDuplicate(qep->V[0],&xi);
551: QEPGetEigenpair(qep,i,&kr,&ki,xr,xi);
552: QEPComputeRelativeError_Private(qep,kr,ki,xr,xi,error);
553: VecDestroy(&xr);
554: VecDestroy(&xi);
555: return(0);
556: }
560: /*@
561: QEPSortEigenvalues - Sorts a list of eigenvalues according to the criterion
562: specified via QEPSetWhichEigenpairs().
564: Not Collective
566: Input Parameters:
567: + qep - the quadratic eigensolver context
568: . n - number of eigenvalues in the list
569: . eigr - pointer to the array containing the eigenvalues
570: - eigi - imaginary part of the eigenvalues (only when using real numbers)
572: Output Parameter:
573: . perm - resulting permutation
575: Note:
576: The result is a list of indices in the original eigenvalue array
577: corresponding to the first nev eigenvalues sorted in the specified
578: criterion.
580: Level: developer
582: .seealso: QEPSortEigenvaluesReal(), QEPSetWhichEigenpairs()
583: @*/
584: PetscErrorCode QEPSortEigenvalues(QEP qep,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm)
585: {
587: PetscScalar re,im;
588: PetscInt i,j,result,tmp;
595: for (i=0; i<n; i++) { perm[i] = i; }
596: /* insertion sort */
597: for (i=n-1; i>=0; i--) {
598: re = eigr[perm[i]];
599: im = eigi[perm[i]];
600: j = i + 1;
601: #if !defined(PETSC_USE_COMPLEX)
602: if (im != 0) {
603: /* complex eigenvalue */
604: i--;
605: im = eigi[perm[i]];
606: }
607: #endif
608: while (j<n) {
609: QEPCompareEigenvalues(qep,re,im,eigr[perm[j]],eigi[perm[j]],&result);
610: if (result >= 0) break;
611: #if !defined(PETSC_USE_COMPLEX)
612: /* keep together every complex conjugated eigenpair */
613: if (im == 0) {
614: if (eigi[perm[j]] == 0) {
615: #endif
616: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = tmp;
617: j++;
618: #if !defined(PETSC_USE_COMPLEX)
619: } else {
620: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp;
621: j+=2;
622: }
623: } else {
624: if (eigi[perm[j]] == 0) {
625: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = perm[j-1]; perm[j-1] = tmp;
626: j++;
627: } else {
628: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = tmp;
629: tmp = perm[j-1]; perm[j-1] = perm[j+1]; perm[j+1] = tmp;
630: j+=2;
631: }
632: }
633: #endif
634: }
635: }
636: return(0);
637: }
641: /*@
642: QEPSortEigenvaluesReal - Sorts a list of eigenvalues according to a certain
643: criterion (version for real eigenvalues only).
645: Not Collective
647: Input Parameters:
648: + qep - the quadratic eigensolver context
649: . n - number of eigenvalue in the list
650: - eig - pointer to the array containing the eigenvalues (real)
652: Output Parameter:
653: . perm - resulting permutation
655: Note:
656: The result is a list of indices in the original eigenvalue array
657: corresponding to the first nev eigenvalues sorted in the specified
658: criterion.
660: Level: developer
662: .seealso: QEPSortEigenvalues(), QEPSetWhichEigenpairs(), QEPCompareEigenvalues()
663: @*/
664: PetscErrorCode QEPSortEigenvaluesReal(QEP qep,PetscInt n,PetscReal *eig,PetscInt *perm)
665: {
667: PetscScalar re;
668: PetscInt i,j,result,tmp;
674: for (i=0; i<n; i++) { perm[i] = i; }
675: /* insertion sort */
676: for (i=1; i<n; i++) {
677: re = eig[perm[i]];
678: j = i-1;
679: QEPCompareEigenvalues(qep,re,0.0,eig[perm[j]],0.0,&result);
680: while (result>0 && j>=0) {
681: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
682: if (j>=0) {
683: QEPCompareEigenvalues(qep,re,0.0,eig[perm[j]],0.0,&result);
684: }
685: }
686: }
687: return(0);
688: }
692: /*@
693: QEPCompareEigenvalues - Compares two (possibly complex) eigenvalues according
694: to a certain criterion.
696: Not Collective
698: Input Parameters:
699: + qep - the quadratic eigensolver context
700: . ar - real part of the 1st eigenvalue
701: . ai - imaginary part of the 1st eigenvalue
702: . br - real part of the 2nd eigenvalue
703: - bi - imaginary part of the 2nd eigenvalue
705: Output Parameter:
706: . res - result of comparison
708: Notes:
709: Returns an integer less than, equal to, or greater than zero if the first
710: eigenvalue is considered to be respectively less than, equal to, or greater
711: than the second one.
713: The criterion of comparison is related to the 'which' parameter set with
714: QEPSetWhichEigenpairs().
716: Level: developer
718: .seealso: QEPSortEigenvalues(), QEPSetWhichEigenpairs()
719: @*/
720: PetscErrorCode QEPCompareEigenvalues(QEP qep,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result)
721: {
722: PetscReal a,b;
727: switch(qep->which) {
728: case QEP_LARGEST_MAGNITUDE:
729: case QEP_SMALLEST_MAGNITUDE:
730: a = SlepcAbsEigenvalue(ar,ai);
731: b = SlepcAbsEigenvalue(br,bi);
732: break;
733: case QEP_LARGEST_REAL:
734: case QEP_SMALLEST_REAL:
735: a = PetscRealPart(ar);
736: b = PetscRealPart(br);
737: break;
738: case QEP_LARGEST_IMAGINARY:
739: case QEP_SMALLEST_IMAGINARY:
740: #if defined(PETSC_USE_COMPLEX)
741: a = PetscImaginaryPart(ar);
742: b = PetscImaginaryPart(br);
743: #else
744: a = PetscAbsReal(ai);
745: b = PetscAbsReal(bi);
746: #endif
747: break;
748: default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of which");
749: }
750: switch(qep->which) {
751: case QEP_LARGEST_MAGNITUDE:
752: case QEP_LARGEST_REAL:
753: case QEP_LARGEST_IMAGINARY:
754: if (a<b) *result = -1;
755: else if (a>b) *result = 1;
756: else *result = 0;
757: break;
758: default:
759: if (a>b) *result = -1;
760: else if (a<b) *result = 1;
761: else *result = 0;
762: }
763: return(0);
764: }
768: /*@
769: QEPGetOperationCounters - Gets the total number of matrix-vector products, dot
770: products, and linear solve iterations used by the QEP object during the last
771: QEPSolve() call.
773: Not Collective
775: Input Parameter:
776: . qep - quadratic eigensolver context
778: Output Parameter:
779: + matvecs - number of matrix-vector product operations
780: . dots - number of dot product operations
781: - lits - number of linear iterations
783: Notes:
784: These counters are reset to zero at each successive call to QEPSolve().
786: Level: intermediate
788: @*/
789: PetscErrorCode QEPGetOperationCounters(QEP qep,PetscInt* matvecs,PetscInt* dots,PetscInt* lits)
790: {
792:
795: if (matvecs) *matvecs = qep->matvecs;
796: if (dots) {
797: if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
798: IPGetOperationCounters(qep->ip,dots);
799: }
800: if (lits) *lits = qep->linits;
801: return(0);
802: }