| 1887 |
jroman |
1 |
/*
|
|
|
2 |
QEP routines related to the solution process.
|
|
|
3 |
|
|
|
4 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
5 |
SLEPc - Scalable Library for Eigenvalue Problem Computations
|
| 2116 |
eromero |
6 |
Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
|
| 1887 |
jroman |
7 |
|
|
|
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.
|
|
|
13 |
|
|
|
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.
|
|
|
18 |
|
|
|
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 |
*/
|
|
|
23 |
|
|
|
24 |
#include "private/qepimpl.h" /*I "slepcqep.h" I*/
|
|
|
25 |
|
|
|
26 |
#undef __FUNCT__
|
|
|
27 |
#define __FUNCT__ "QEPSolve"
|
|
|
28 |
/*@
|
|
|
29 |
QEPSolve - Solves the quadratic eigensystem.
|
|
|
30 |
|
|
|
31 |
Collective on QEP
|
|
|
32 |
|
|
|
33 |
Input Parameter:
|
|
|
34 |
. qep - eigensolver context obtained from QEPCreate()
|
|
|
35 |
|
|
|
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
|
|
|
40 |
|
|
|
41 |
Level: beginner
|
|
|
42 |
|
|
|
43 |
.seealso: QEPCreate(), QEPSetUp(), QEPDestroy(), QEPSetTolerances()
|
|
|
44 |
@*/
|
|
|
45 |
PetscErrorCode QEPSolve(QEP qep)
|
|
|
46 |
{
|
|
|
47 |
PetscErrorCode ierr;
|
|
|
48 |
PetscInt i;
|
|
|
49 |
PetscReal re,im;
|
|
|
50 |
PetscTruth flg;
|
|
|
51 |
PetscViewer viewer;
|
|
|
52 |
PetscDraw draw;
|
|
|
53 |
PetscDrawSP drawsp;
|
|
|
54 |
char filename[PETSC_MAX_PATH_LEN];
|
|
|
55 |
|
|
|
56 |
PetscFunctionBegin;
|
| 2213 |
jroman |
57 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
58 |
|
|
|
59 |
flg = PETSC_FALSE;
|
|
|
60 |
ierr = PetscOptionsGetTruth(((PetscObject)qep)->prefix,"-qep_view_binary",&flg,PETSC_NULL);CHKERRQ(ierr);
|
|
|
61 |
if (flg) {
|
|
|
62 |
ierr = MatView(qep->M,PETSC_VIEWER_BINARY_(((PetscObject)qep)->comm));CHKERRQ(ierr);
|
|
|
63 |
ierr = MatView(qep->C,PETSC_VIEWER_BINARY_(((PetscObject)qep)->comm));CHKERRQ(ierr);
|
|
|
64 |
ierr = MatView(qep->K,PETSC_VIEWER_BINARY_(((PetscObject)qep)->comm));CHKERRQ(ierr);
|
|
|
65 |
}
|
|
|
66 |
|
|
|
67 |
/* reset the convergence flag from the previous solves */
|
|
|
68 |
qep->reason = QEP_CONVERGED_ITERATING;
|
|
|
69 |
|
|
|
70 |
if (!qep->setupcalled){ ierr = QEPSetUp(qep);CHKERRQ(ierr); }
|
|
|
71 |
ierr = IPResetOperationCounters(qep->ip);CHKERRQ(ierr);
|
|
|
72 |
qep->nconv = 0;
|
|
|
73 |
qep->its = 0;
|
| 1896 |
jroman |
74 |
qep->matvecs = 0;
|
|
|
75 |
qep->linits = 0;
|
| 1887 |
jroman |
76 |
for (i=0;i<qep->ncv;i++) qep->eigr[i]=qep->eigi[i]=qep->errest[i]=0.0;
|
|
|
77 |
QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,qep->ncv);
|
|
|
78 |
|
|
|
79 |
ierr = PetscLogEventBegin(QEP_Solve,qep,0,0,0);CHKERRQ(ierr);
|
|
|
80 |
ierr = (*qep->ops->solve)(qep);CHKERRQ(ierr);
|
|
|
81 |
ierr = PetscLogEventEnd(QEP_Solve,qep,0,0,0);CHKERRQ(ierr);
|
|
|
82 |
|
|
|
83 |
if (!qep->reason) {
|
|
|
84 |
SETERRQ(1,"Internal error, solver returned without setting converged reason");
|
|
|
85 |
}
|
|
|
86 |
|
|
|
87 |
#ifndef PETSC_USE_COMPLEX
|
|
|
88 |
/* reorder conjugate eigenvalues (positive imaginary first) */
|
|
|
89 |
for (i=0;i<qep->nconv-1;i++) {
|
|
|
90 |
if (qep->eigi[i] != 0) {
|
|
|
91 |
if (qep->eigi[i] < 0) {
|
|
|
92 |
qep->eigi[i] = -qep->eigi[i];
|
|
|
93 |
qep->eigi[i+1] = -qep->eigi[i+1];
|
|
|
94 |
ierr = VecScale(qep->V[i+1],-1.0);CHKERRQ(ierr);
|
|
|
95 |
}
|
|
|
96 |
i++;
|
|
|
97 |
}
|
|
|
98 |
}
|
|
|
99 |
#endif
|
|
|
100 |
|
|
|
101 |
/* sort eigenvalues according to qep->which parameter */
|
|
|
102 |
ierr = PetscFree(qep->perm);CHKERRQ(ierr);
|
|
|
103 |
if (qep->nconv > 0) {
|
|
|
104 |
ierr = PetscMalloc(sizeof(PetscInt)*qep->nconv,&qep->perm);CHKERRQ(ierr);
|
|
|
105 |
ierr = QEPSortEigenvalues(qep,qep->nconv,qep->eigr,qep->eigi,qep->perm);CHKERRQ(ierr);
|
|
|
106 |
}
|
|
|
107 |
|
|
|
108 |
ierr = PetscOptionsGetString(((PetscObject)qep)->prefix,"-qep_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
|
|
|
109 |
if (flg && !PetscPreLoadingOn) {
|
|
|
110 |
ierr = PetscViewerASCIIOpen(((PetscObject)qep)->comm,filename,&viewer);CHKERRQ(ierr);
|
|
|
111 |
ierr = QEPView(qep,viewer);CHKERRQ(ierr);
|
|
|
112 |
ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
|
|
|
113 |
}
|
|
|
114 |
|
|
|
115 |
flg = PETSC_FALSE;
|
|
|
116 |
ierr = PetscOptionsGetTruth(((PetscObject)qep)->prefix,"-qep_plot_eigs",&flg,PETSC_NULL);CHKERRQ(ierr);
|
|
|
117 |
if (flg) {
|
|
|
118 |
ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
|
|
|
119 |
PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);CHKERRQ(ierr);
|
|
|
120 |
ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
|
|
|
121 |
ierr = PetscDrawSPCreate(draw,1,&drawsp);CHKERRQ(ierr);
|
|
|
122 |
for(i=0;i<qep->nconv;i++) {
|
|
|
123 |
#if defined(PETSC_USE_COMPLEX)
|
|
|
124 |
re = PetscRealPart(qep->eigr[i]);
|
|
|
125 |
im = PetscImaginaryPart(qep->eigi[i]);
|
|
|
126 |
#else
|
|
|
127 |
re = qep->eigr[i];
|
|
|
128 |
im = qep->eigi[i];
|
|
|
129 |
#endif
|
|
|
130 |
ierr = PetscDrawSPAddPoint(drawsp,&re,&im);CHKERRQ(ierr);
|
|
|
131 |
}
|
|
|
132 |
ierr = PetscDrawSPDraw(drawsp);CHKERRQ(ierr);
|
|
|
133 |
ierr = PetscDrawSPDestroy(drawsp);CHKERRQ(ierr);
|
|
|
134 |
ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
|
|
|
135 |
}
|
|
|
136 |
|
| 2080 |
eromero |
137 |
/* Remove the initial subspace */
|
|
|
138 |
qep->nini = 0;
|
|
|
139 |
|
| 1887 |
jroman |
140 |
PetscFunctionReturn(0);
|
|
|
141 |
}
|
|
|
142 |
|
|
|
143 |
#undef __FUNCT__
|
|
|
144 |
#define __FUNCT__ "QEPGetIterationNumber"
|
|
|
145 |
/*@
|
|
|
146 |
QEPGetIterationNumber - Gets the current iteration number. If the
|
|
|
147 |
call to QEPSolve() is complete, then it returns the number of iterations
|
|
|
148 |
carried out by the solution method.
|
|
|
149 |
|
|
|
150 |
Not Collective
|
|
|
151 |
|
|
|
152 |
Input Parameter:
|
|
|
153 |
. qep - the quadratic eigensolver context
|
|
|
154 |
|
|
|
155 |
Output Parameter:
|
|
|
156 |
. its - number of iterations
|
|
|
157 |
|
|
|
158 |
Level: intermediate
|
|
|
159 |
|
|
|
160 |
Note:
|
|
|
161 |
During the i-th iteration this call returns i-1. If QEPSolve() is
|
|
|
162 |
complete, then parameter "its" contains either the iteration number at
|
|
|
163 |
which convergence was successfully reached, or failure was detected.
|
|
|
164 |
Call QEPGetConvergedReason() to determine if the solver converged or
|
|
|
165 |
failed and why.
|
|
|
166 |
|
|
|
167 |
.seealso: QEPGetConvergedReason(), QEPSetTolerances()
|
|
|
168 |
@*/
|
|
|
169 |
PetscErrorCode QEPGetIterationNumber(QEP qep,PetscInt *its)
|
|
|
170 |
{
|
|
|
171 |
PetscFunctionBegin;
|
| 2213 |
jroman |
172 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
173 |
PetscValidIntPointer(its,2);
|
|
|
174 |
*its = qep->its;
|
|
|
175 |
PetscFunctionReturn(0);
|
|
|
176 |
}
|
|
|
177 |
|
|
|
178 |
#undef __FUNCT__
|
|
|
179 |
#define __FUNCT__ "QEPGetConverged"
|
|
|
180 |
/*@
|
|
|
181 |
QEPGetConverged - Gets the number of converged eigenpairs.
|
|
|
182 |
|
|
|
183 |
Not Collective
|
|
|
184 |
|
|
|
185 |
Input Parameter:
|
|
|
186 |
. qep - the quadratic eigensolver context
|
|
|
187 |
|
|
|
188 |
Output Parameter:
|
|
|
189 |
. nconv - number of converged eigenpairs
|
|
|
190 |
|
|
|
191 |
Note:
|
|
|
192 |
This function should be called after QEPSolve() has finished.
|
|
|
193 |
|
|
|
194 |
Level: beginner
|
|
|
195 |
|
|
|
196 |
.seealso: QEPSetDimensions(), QEPSolve()
|
|
|
197 |
@*/
|
|
|
198 |
PetscErrorCode QEPGetConverged(QEP qep,PetscInt *nconv)
|
|
|
199 |
{
|
|
|
200 |
PetscFunctionBegin;
|
| 2213 |
jroman |
201 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
202 |
PetscValidIntPointer(nconv,2);
|
|
|
203 |
*nconv = qep->nconv;
|
|
|
204 |
PetscFunctionReturn(0);
|
|
|
205 |
}
|
|
|
206 |
|
|
|
207 |
|
|
|
208 |
#undef __FUNCT__
|
|
|
209 |
#define __FUNCT__ "QEPGetConvergedReason"
|
|
|
210 |
/*@C
|
|
|
211 |
QEPGetConvergedReason - Gets the reason why the QEPSolve() iteration was
|
|
|
212 |
stopped.
|
|
|
213 |
|
|
|
214 |
Not Collective
|
|
|
215 |
|
|
|
216 |
Input Parameter:
|
|
|
217 |
. qep - the quadratic eigensolver context
|
|
|
218 |
|
|
|
219 |
Output Parameter:
|
|
|
220 |
. reason - negative value indicates diverged, positive value converged
|
|
|
221 |
|
|
|
222 |
Possible values for reason:
|
|
|
223 |
+ QEP_CONVERGED_TOL - converged up to tolerance
|
|
|
224 |
. QEP_DIVERGED_ITS - required more than its to reach convergence
|
|
|
225 |
- QEP_DIVERGED_BREAKDOWN - generic breakdown in method
|
|
|
226 |
|
|
|
227 |
Note:
|
|
|
228 |
Can only be called after the call to QEPSolve() is complete.
|
|
|
229 |
|
|
|
230 |
Level: intermediate
|
|
|
231 |
|
|
|
232 |
.seealso: QEPSetTolerances(), QEPSolve(), QEPConvergedReason
|
|
|
233 |
@*/
|
|
|
234 |
PetscErrorCode QEPGetConvergedReason(QEP qep,QEPConvergedReason *reason)
|
|
|
235 |
{
|
|
|
236 |
PetscFunctionBegin;
|
| 2213 |
jroman |
237 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
238 |
PetscValidIntPointer(reason,2);
|
|
|
239 |
*reason = qep->reason;
|
|
|
240 |
PetscFunctionReturn(0);
|
|
|
241 |
}
|
|
|
242 |
|
|
|
243 |
#undef __FUNCT__
|
|
|
244 |
#define __FUNCT__ "QEPGetEigenpair"
|
|
|
245 |
/*@
|
|
|
246 |
QEPGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
|
|
|
247 |
QEPSolve(). The solution consists in both the eigenvalue and the eigenvector.
|
|
|
248 |
|
|
|
249 |
Not Collective, but vectors are shared by all processors that share the QEP
|
|
|
250 |
|
|
|
251 |
Input Parameters:
|
|
|
252 |
+ qep - quadratic eigensolver context
|
|
|
253 |
- i - index of the solution
|
|
|
254 |
|
|
|
255 |
Output Parameters:
|
|
|
256 |
+ eigr - real part of eigenvalue
|
|
|
257 |
. eigi - imaginary part of eigenvalue
|
|
|
258 |
. Vr - real part of eigenvector
|
|
|
259 |
- Vi - imaginary part of eigenvector
|
|
|
260 |
|
|
|
261 |
Notes:
|
|
|
262 |
If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
|
|
|
263 |
configured with complex scalars the eigenvalue is stored
|
|
|
264 |
directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
|
|
|
265 |
set to zero).
|
|
|
266 |
|
|
|
267 |
The index i should be a value between 0 and nconv-1 (see QEPGetConverged()).
|
|
|
268 |
Eigenpairs are indexed according to the ordering criterion established
|
|
|
269 |
with QEPSetWhichEigenpairs().
|
|
|
270 |
|
|
|
271 |
Level: beginner
|
|
|
272 |
|
|
|
273 |
.seealso: QEPSolve(), QEPGetConverged(), QEPSetWhichEigenpairs()
|
|
|
274 |
@*/
|
|
|
275 |
PetscErrorCode QEPGetEigenpair(QEP qep, PetscInt i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
|
|
|
276 |
{
|
|
|
277 |
PetscInt k;
|
|
|
278 |
PetscErrorCode ierr;
|
|
|
279 |
|
|
|
280 |
PetscFunctionBegin;
|
| 2213 |
jroman |
281 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
282 |
if (!qep->eigr || !qep->eigi || !qep->V) {
|
|
|
283 |
SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "QEPSolve must be called first");
|
|
|
284 |
}
|
|
|
285 |
if (i<0 || i>=qep->nconv) {
|
|
|
286 |
SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
|
|
|
287 |
}
|
|
|
288 |
|
|
|
289 |
if (!qep->perm) k = i;
|
|
|
290 |
else k = qep->perm[i];
|
|
|
291 |
|
|
|
292 |
/* eigenvalue */
|
|
|
293 |
#ifdef PETSC_USE_COMPLEX
|
|
|
294 |
if (eigr) *eigr = qep->eigr[k];
|
|
|
295 |
if (eigi) *eigi = 0;
|
|
|
296 |
#else
|
|
|
297 |
if (eigr) *eigr = qep->eigr[k];
|
|
|
298 |
if (eigi) *eigi = qep->eigi[k];
|
|
|
299 |
#endif
|
|
|
300 |
|
|
|
301 |
/* eigenvector */
|
|
|
302 |
#ifdef PETSC_USE_COMPLEX
|
|
|
303 |
if (Vr) { ierr = VecCopy(qep->V[k],Vr);CHKERRQ(ierr); }
|
|
|
304 |
if (Vi) { ierr = VecSet(Vi,0.0);CHKERRQ(ierr); }
|
|
|
305 |
#else
|
|
|
306 |
if (qep->eigi[k]>0) { /* first value of conjugate pair */
|
|
|
307 |
if (Vr) { ierr = VecCopy(qep->V[k],Vr);CHKERRQ(ierr); }
|
|
|
308 |
if (Vi) { ierr = VecCopy(qep->V[k+1],Vi);CHKERRQ(ierr); }
|
|
|
309 |
} else if (qep->eigi[k]<0) { /* second value of conjugate pair */
|
|
|
310 |
if (Vr) { ierr = VecCopy(qep->V[k-1],Vr);CHKERRQ(ierr); }
|
|
|
311 |
if (Vi) {
|
|
|
312 |
ierr = VecCopy(qep->V[k],Vi);CHKERRQ(ierr);
|
|
|
313 |
ierr = VecScale(Vi,-1.0);CHKERRQ(ierr);
|
|
|
314 |
}
|
|
|
315 |
} else { /* real eigenvalue */
|
|
|
316 |
if (Vr) { ierr = VecCopy(qep->V[k],Vr);CHKERRQ(ierr); }
|
|
|
317 |
if (Vi) { ierr = VecSet(Vi,0.0);CHKERRQ(ierr); }
|
|
|
318 |
}
|
|
|
319 |
#endif
|
|
|
320 |
|
|
|
321 |
PetscFunctionReturn(0);
|
|
|
322 |
}
|
|
|
323 |
|
|
|
324 |
#undef __FUNCT__
|
|
|
325 |
#define __FUNCT__ "QEPGetErrorEstimate"
|
|
|
326 |
/*@
|
|
|
327 |
QEPGetErrorEstimate - Returns the error estimate associated to the i-th
|
|
|
328 |
computed eigenpair.
|
|
|
329 |
|
|
|
330 |
Not Collective
|
|
|
331 |
|
|
|
332 |
Input Parameter:
|
|
|
333 |
+ qep - quadratic eigensolver context
|
|
|
334 |
- i - index of eigenpair
|
|
|
335 |
|
|
|
336 |
Output Parameter:
|
|
|
337 |
. errest - the error estimate
|
|
|
338 |
|
|
|
339 |
Notes:
|
|
|
340 |
This is the error estimate used internally by the eigensolver. The actual
|
|
|
341 |
error bound can be computed with QEPComputeRelativeError(). See also the users
|
|
|
342 |
manual for details.
|
|
|
343 |
|
|
|
344 |
Level: advanced
|
|
|
345 |
|
|
|
346 |
.seealso: QEPComputeRelativeError()
|
|
|
347 |
@*/
|
|
|
348 |
PetscErrorCode QEPGetErrorEstimate(QEP qep, PetscInt i, PetscReal *errest)
|
|
|
349 |
{
|
|
|
350 |
PetscFunctionBegin;
|
| 2213 |
jroman |
351 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
352 |
if (!qep->eigr || !qep->eigi) {
|
|
|
353 |
SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "QEPSolve must be called first");
|
|
|
354 |
}
|
|
|
355 |
if (i<0 || i>=qep->nconv) {
|
|
|
356 |
SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
|
|
|
357 |
}
|
|
|
358 |
if (qep->perm) i = qep->perm[i];
|
|
|
359 |
if (errest) *errest = qep->errest[i];
|
|
|
360 |
PetscFunctionReturn(0);
|
|
|
361 |
}
|
|
|
362 |
|
|
|
363 |
#undef __FUNCT__
|
|
|
364 |
#define __FUNCT__ "QEPComputeResidualNorm_Private"
|
|
|
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 |
{
|
|
|
371 |
PetscErrorCode ierr;
|
|
|
372 |
Vec u,w;
|
| 1888 |
jroman |
373 |
Mat M=qep->M,C=qep->C,K=qep->K;
|
| 1887 |
jroman |
374 |
#ifndef PETSC_USE_COMPLEX
|
|
|
375 |
Vec v,y,z;
|
|
|
376 |
PetscReal ni,nr;
|
|
|
377 |
PetscScalar a1,a2;
|
|
|
378 |
#endif
|
|
|
379 |
|
|
|
380 |
PetscFunctionBegin;
|
| 1952 |
jroman |
381 |
ierr = VecDuplicate(qep->V[0],&u);CHKERRQ(ierr);
|
| 1887 |
jroman |
382 |
ierr = VecDuplicate(u,&w);CHKERRQ(ierr);
|
|
|
383 |
|
|
|
384 |
#ifndef PETSC_USE_COMPLEX
|
|
|
385 |
if (ki == 0 ||
|
|
|
386 |
PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
|
|
|
387 |
#endif
|
| 1888 |
jroman |
388 |
ierr = MatMult(K,xr,u);CHKERRQ(ierr); /* u=K*x */
|
| 1887 |
jroman |
389 |
if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
|
|
|
390 |
ierr = MatMult(C,xr,w);CHKERRQ(ierr); /* w=C*x */
|
| 1888 |
jroman |
391 |
ierr = VecAXPY(u,kr,w);CHKERRQ(ierr); /* u=l*C*x+K*x */
|
| 1893 |
jroman |
392 |
ierr = MatMult(M,xr,w);CHKERRQ(ierr); /* w=M*x */
|
| 1888 |
jroman |
393 |
ierr = VecAXPY(u,kr*kr,w);CHKERRQ(ierr); /* u=l^2*M*x+l*C*x+K*x */
|
| 1887 |
jroman |
394 |
}
|
|
|
395 |
ierr = VecNorm(u,NORM_2,norm);CHKERRQ(ierr);
|
|
|
396 |
#ifndef PETSC_USE_COMPLEX
|
|
|
397 |
} else {
|
|
|
398 |
ierr = VecDuplicate(u,&v);CHKERRQ(ierr);
|
|
|
399 |
ierr = VecDuplicate(u,&y);CHKERRQ(ierr);
|
|
|
400 |
ierr = VecDuplicate(u,&z);CHKERRQ(ierr);
|
|
|
401 |
a1 = kr*kr-ki*ki;
|
|
|
402 |
a2 = 2.0*kr*ki;
|
| 1893 |
jroman |
403 |
ierr = MatMult(K,xr,u);CHKERRQ(ierr); /* u=K*xr */
|
|
|
404 |
if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
|
|
|
405 |
ierr = MatMult(C,xr,v);CHKERRQ(ierr); /* v=C*xr */
|
|
|
406 |
ierr = MatMult(C,xi,w);CHKERRQ(ierr); /* w=C*xi */
|
|
|
407 |
ierr = MatMult(M,xr,y);CHKERRQ(ierr); /* y=M*xr */
|
|
|
408 |
ierr = MatMult(M,xi,z);CHKERRQ(ierr); /* z=M*xi */
|
|
|
409 |
ierr = VecAXPY(u,kr,v);CHKERRQ(ierr); /* u=kr*C*xr+K*xr */
|
|
|
410 |
ierr = VecAXPY(u,-ki,w);CHKERRQ(ierr); /* u=kr*C*xr-ki*C*xi+K*xr */
|
|
|
411 |
ierr = VecAXPY(u,a1,y);CHKERRQ(ierr); /* u=a1*M*xr+kr*C*xr-ki*C*xi+K*xr */
|
|
|
412 |
ierr = VecAXPY(u,-a2,z);CHKERRQ(ierr); /* u=a1*M*xr-a2*M*xi+kr*C*xr-ki*C*xi+K*xr */
|
|
|
413 |
}
|
| 1887 |
jroman |
414 |
ierr = VecNorm(u,NORM_2,&nr);CHKERRQ(ierr);
|
| 1888 |
jroman |
415 |
ierr = MatMult(K,xi,u);CHKERRQ(ierr); /* u=K*xi */
|
| 1893 |
jroman |
416 |
if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
|
|
|
417 |
ierr = VecAXPY(u,kr,w);CHKERRQ(ierr); /* u=kr*C*xi+K*xi */
|
|
|
418 |
ierr = VecAXPY(u,ki,v);CHKERRQ(ierr); /* u=kr*C*xi+ki*C*xi+K*xi */
|
|
|
419 |
ierr = VecAXPY(u,a1,z);CHKERRQ(ierr); /* u=a1*M*xi+kr*C*xi+ki*C*xi+K*xi */
|
|
|
420 |
ierr = VecAXPY(u,a2,y);CHKERRQ(ierr); /* u=a1*M*xi+a2*M*ki+kr*C*xi+ki*C*xi+K*xi */
|
|
|
421 |
}
|
| 1887 |
jroman |
422 |
ierr = VecNorm(u,NORM_2,&ni);CHKERRQ(ierr);
|
|
|
423 |
*norm = SlepcAbsEigenvalue(nr,ni);
|
|
|
424 |
ierr = VecDestroy(v);CHKERRQ(ierr);
|
|
|
425 |
ierr = VecDestroy(y);CHKERRQ(ierr);
|
|
|
426 |
ierr = VecDestroy(z);CHKERRQ(ierr);
|
|
|
427 |
}
|
|
|
428 |
#endif
|
|
|
429 |
|
|
|
430 |
ierr = VecDestroy(w);CHKERRQ(ierr);
|
|
|
431 |
ierr = VecDestroy(u);CHKERRQ(ierr);
|
|
|
432 |
PetscFunctionReturn(0);
|
|
|
433 |
}
|
|
|
434 |
|
|
|
435 |
#undef __FUNCT__
|
|
|
436 |
#define __FUNCT__ "QEPComputeResidualNorm"
|
|
|
437 |
/*@
|
|
|
438 |
QEPComputeResidualNorm - Computes the norm of the residual vector associated with
|
|
|
439 |
the i-th computed eigenpair.
|
|
|
440 |
|
|
|
441 |
Collective on QEP
|
|
|
442 |
|
|
|
443 |
Input Parameter:
|
|
|
444 |
. qep - the quadratic eigensolver context
|
|
|
445 |
. i - the solution index
|
|
|
446 |
|
|
|
447 |
Output Parameter:
|
| 1888 |
jroman |
448 |
. norm - the residual norm, computed as ||(l^2*M+l*C+K)x||_2 where l is the
|
| 1887 |
jroman |
449 |
eigenvalue and x is the eigenvector.
|
| 1893 |
jroman |
450 |
If l=0 then the residual norm is computed as ||Kx||_2.
|
| 1887 |
jroman |
451 |
|
|
|
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().
|
|
|
456 |
|
|
|
457 |
Level: beginner
|
|
|
458 |
|
|
|
459 |
.seealso: QEPSolve(), QEPGetConverged(), QEPSetWhichEigenpairs()
|
|
|
460 |
@*/
|
|
|
461 |
PetscErrorCode QEPComputeResidualNorm(QEP qep, PetscInt i, PetscReal *norm)
|
|
|
462 |
{
|
|
|
463 |
PetscErrorCode ierr;
|
|
|
464 |
Vec xr,xi;
|
|
|
465 |
PetscScalar kr,ki;
|
|
|
466 |
|
|
|
467 |
PetscFunctionBegin;
|
| 2213 |
jroman |
468 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
469 |
PetscValidPointer(norm,3);
|
| 1952 |
jroman |
470 |
ierr = VecDuplicate(qep->V[0],&xr);CHKERRQ(ierr);
|
|
|
471 |
ierr = VecDuplicate(qep->V[0],&xi);CHKERRQ(ierr);
|
| 1887 |
jroman |
472 |
ierr = QEPGetEigenpair(qep,i,&kr,&ki,xr,xi);CHKERRQ(ierr);
|
|
|
473 |
ierr = QEPComputeResidualNorm_Private(qep,kr,ki,xr,xi,norm);CHKERRQ(ierr);
|
|
|
474 |
ierr = VecDestroy(xr);CHKERRQ(ierr);
|
|
|
475 |
ierr = VecDestroy(xi);CHKERRQ(ierr);
|
|
|
476 |
PetscFunctionReturn(0);
|
|
|
477 |
}
|
|
|
478 |
|
|
|
479 |
#undef __FUNCT__
|
|
|
480 |
#define __FUNCT__ "QEPComputeRelativeError_Private"
|
|
|
481 |
/*
|
|
|
482 |
QEPComputeRelativeError_Private - Computes the relative error bound
|
|
|
483 |
associated with an eigenpair.
|
|
|
484 |
*/
|
|
|
485 |
PetscErrorCode QEPComputeRelativeError_Private(QEP qep, PetscScalar kr, PetscScalar ki, Vec xr, Vec xi, PetscReal *error)
|
|
|
486 |
{
|
|
|
487 |
PetscErrorCode ierr;
|
|
|
488 |
PetscReal norm, er;
|
|
|
489 |
#ifndef PETSC_USE_COMPLEX
|
|
|
490 |
PetscReal ei;
|
|
|
491 |
#endif
|
|
|
492 |
|
|
|
493 |
PetscFunctionBegin;
|
| 1893 |
jroman |
494 |
ierr = QEPComputeResidualNorm_Private(qep,kr,ki,xr,xi,&norm);CHKERRQ(ierr);
|
| 1887 |
jroman |
495 |
|
|
|
496 |
#ifndef PETSC_USE_COMPLEX
|
|
|
497 |
if (ki == 0 ||
|
|
|
498 |
PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
|
|
|
499 |
#endif
|
| 1893 |
jroman |
500 |
ierr = VecNorm(xr,NORM_2,&er);CHKERRQ(ierr);
|
| 1887 |
jroman |
501 |
if (PetscAbsScalar(kr) > norm) {
|
| 1893 |
jroman |
502 |
*error = norm/(PetscAbsScalar(kr)*er);
|
| 1887 |
jroman |
503 |
} else {
|
| 1893 |
jroman |
504 |
*error = norm/er;
|
| 1887 |
jroman |
505 |
}
|
|
|
506 |
#ifndef PETSC_USE_COMPLEX
|
|
|
507 |
} else {
|
| 1893 |
jroman |
508 |
ierr = VecNorm(xr,NORM_2,&er);CHKERRQ(ierr);
|
|
|
509 |
ierr = VecNorm(xi,NORM_2,&ei);CHKERRQ(ierr);
|
| 1887 |
jroman |
510 |
if (SlepcAbsEigenvalue(kr,ki) > norm) {
|
| 1893 |
jroman |
511 |
*error = norm/(SlepcAbsEigenvalue(kr,ki)*SlepcAbsEigenvalue(er,ei));
|
| 1887 |
jroman |
512 |
} else {
|
| 1893 |
jroman |
513 |
*error = norm/SlepcAbsEigenvalue(er,ei);
|
| 1887 |
jroman |
514 |
}
|
|
|
515 |
}
|
|
|
516 |
#endif
|
|
|
517 |
|
|
|
518 |
PetscFunctionReturn(0);
|
|
|
519 |
}
|
|
|
520 |
|
|
|
521 |
#undef __FUNCT__
|
|
|
522 |
#define __FUNCT__ "QEPComputeRelativeError"
|
|
|
523 |
/*@
|
|
|
524 |
QEPComputeRelativeError - Computes the relative error bound associated
|
|
|
525 |
with the i-th computed eigenpair.
|
|
|
526 |
|
|
|
527 |
Collective on QEP
|
|
|
528 |
|
|
|
529 |
Input Parameter:
|
|
|
530 |
. qep - the quadratic eigensolver context
|
|
|
531 |
. i - the solution index
|
|
|
532 |
|
|
|
533 |
Output Parameter:
|
| 1888 |
jroman |
534 |
. error - the relative error bound, computed as ||(l^2*M+l*C+K)x||_2/||lx||_2 where
|
| 1887 |
jroman |
535 |
l is the eigenvalue and x is the eigenvector.
|
| 1893 |
jroman |
536 |
If l=0 the relative error is computed as ||Kx||_2/||x||_2.
|
| 1887 |
jroman |
537 |
|
|
|
538 |
Level: beginner
|
|
|
539 |
|
|
|
540 |
.seealso: QEPSolve(), QEPComputeResidualNorm(), QEPGetErrorEstimate()
|
|
|
541 |
@*/
|
|
|
542 |
PetscErrorCode QEPComputeRelativeError(QEP qep, PetscInt i, PetscReal *error)
|
|
|
543 |
{
|
|
|
544 |
PetscErrorCode ierr;
|
|
|
545 |
Vec xr,xi;
|
|
|
546 |
PetscScalar kr,ki;
|
|
|
547 |
|
|
|
548 |
PetscFunctionBegin;
|
| 2213 |
jroman |
549 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
550 |
PetscValidPointer(error,3);
|
| 1952 |
jroman |
551 |
ierr = VecDuplicate(qep->V[0],&xr);CHKERRQ(ierr);
|
|
|
552 |
ierr = VecDuplicate(qep->V[0],&xi);CHKERRQ(ierr);
|
| 1887 |
jroman |
553 |
ierr = QEPGetEigenpair(qep,i,&kr,&ki,xr,xi);CHKERRQ(ierr);
|
|
|
554 |
ierr = QEPComputeRelativeError_Private(qep,kr,ki,xr,xi,error);CHKERRQ(ierr);
|
|
|
555 |
ierr = VecDestroy(xr);CHKERRQ(ierr);
|
|
|
556 |
ierr = VecDestroy(xi);CHKERRQ(ierr);
|
|
|
557 |
PetscFunctionReturn(0);
|
|
|
558 |
}
|
|
|
559 |
|
|
|
560 |
#undef __FUNCT__
|
|
|
561 |
#define __FUNCT__ "QEPSortEigenvalues"
|
|
|
562 |
/*@
|
|
|
563 |
QEPSortEigenvalues - Sorts a list of eigenvalues according to the criterion
|
|
|
564 |
specified via QEPSetWhichEigenpairs().
|
|
|
565 |
|
|
|
566 |
Not Collective
|
|
|
567 |
|
|
|
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)
|
|
|
573 |
|
|
|
574 |
Output Parameter:
|
|
|
575 |
. perm - resulting permutation
|
|
|
576 |
|
|
|
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.
|
|
|
581 |
|
|
|
582 |
Level: developer
|
|
|
583 |
|
|
|
584 |
.seealso: QEPSortEigenvaluesReal(), QEPSetWhichEigenpairs()
|
|
|
585 |
@*/
|
|
|
586 |
PetscErrorCode QEPSortEigenvalues(QEP qep,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm)
|
|
|
587 |
{
|
|
|
588 |
PetscErrorCode ierr;
|
|
|
589 |
PetscScalar re,im;
|
|
|
590 |
PetscInt i,j,result,tmp;
|
|
|
591 |
|
|
|
592 |
PetscFunctionBegin;
|
|
|
593 |
for (i=0; i<n; i++) { perm[i] = i; }
|
|
|
594 |
/* insertion sort */
|
|
|
595 |
for (i=n-1; i>=0; i--) {
|
|
|
596 |
re = eigr[perm[i]];
|
|
|
597 |
im = eigi[perm[i]];
|
|
|
598 |
j = i + 1;
|
|
|
599 |
#ifndef PETSC_USE_COMPLEX
|
|
|
600 |
if (im != 0) {
|
|
|
601 |
/* complex eigenvalue */
|
|
|
602 |
i--;
|
|
|
603 |
im = eigi[perm[i]];
|
|
|
604 |
}
|
|
|
605 |
#endif
|
|
|
606 |
while (j<n) {
|
|
|
607 |
ierr = QEPCompareEigenvalues(qep,re,im,eigr[perm[j]],eigi[perm[j]],&result);CHKERRQ(ierr);
|
|
|
608 |
if (result >= 0) break;
|
|
|
609 |
#ifndef PETSC_USE_COMPLEX
|
|
|
610 |
/* keep together every complex conjugated eigenpair */
|
|
|
611 |
if (im == 0) {
|
|
|
612 |
if (eigi[perm[j]] == 0) {
|
|
|
613 |
#endif
|
|
|
614 |
tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = tmp;
|
|
|
615 |
j++;
|
|
|
616 |
#ifndef PETSC_USE_COMPLEX
|
|
|
617 |
} else {
|
|
|
618 |
tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp;
|
|
|
619 |
j+=2;
|
|
|
620 |
}
|
|
|
621 |
} else {
|
|
|
622 |
if (eigi[perm[j]] == 0) {
|
|
|
623 |
tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = perm[j-1]; perm[j-1] = tmp;
|
|
|
624 |
j++;
|
|
|
625 |
} else {
|
|
|
626 |
tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = tmp;
|
|
|
627 |
tmp = perm[j-1]; perm[j-1] = perm[j+1]; perm[j+1] = tmp;
|
|
|
628 |
j+=2;
|
|
|
629 |
}
|
|
|
630 |
}
|
|
|
631 |
#endif
|
|
|
632 |
}
|
|
|
633 |
}
|
|
|
634 |
PetscFunctionReturn(0);
|
|
|
635 |
}
|
|
|
636 |
|
|
|
637 |
#undef __FUNCT__
|
|
|
638 |
#define __FUNCT__ "QEPSortEigenvaluesReal"
|
|
|
639 |
/*@
|
|
|
640 |
QEPSortEigenvaluesReal - Sorts a list of eigenvalues according to a certain
|
|
|
641 |
criterion (version for real eigenvalues only).
|
|
|
642 |
|
|
|
643 |
Not Collective
|
|
|
644 |
|
|
|
645 |
Input Parameters:
|
|
|
646 |
+ qep - the quadratic eigensolver context
|
|
|
647 |
. n - number of eigenvalue in the list
|
|
|
648 |
- eig - pointer to the array containing the eigenvalues (real)
|
|
|
649 |
|
|
|
650 |
Output Parameter:
|
|
|
651 |
. perm - resulting permutation
|
|
|
652 |
|
|
|
653 |
Note:
|
|
|
654 |
The result is a list of indices in the original eigenvalue array
|
|
|
655 |
corresponding to the first nev eigenvalues sorted in the specified
|
|
|
656 |
criterion.
|
|
|
657 |
|
|
|
658 |
Level: developer
|
|
|
659 |
|
|
|
660 |
.seealso: QEPSortEigenvalues(), QEPSetWhichEigenpairs(), QEPCompareEigenvalues()
|
|
|
661 |
@*/
|
|
|
662 |
PetscErrorCode QEPSortEigenvaluesReal(QEP qep,PetscInt n,PetscReal *eig,PetscInt *perm)
|
|
|
663 |
{
|
|
|
664 |
PetscErrorCode ierr;
|
|
|
665 |
PetscScalar re;
|
|
|
666 |
PetscInt i,j,result,tmp;
|
|
|
667 |
|
|
|
668 |
PetscFunctionBegin;
|
|
|
669 |
for (i=0; i<n; i++) { perm[i] = i; }
|
|
|
670 |
/* insertion sort */
|
|
|
671 |
for (i=1; i<n; i++) {
|
|
|
672 |
re = eig[perm[i]];
|
|
|
673 |
j = i-1;
|
|
|
674 |
ierr = QEPCompareEigenvalues(qep,re,0.0,eig[perm[j]],0.0,&result);CHKERRQ(ierr);
|
|
|
675 |
while (result>0 && j>=0) {
|
|
|
676 |
tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
|
|
|
677 |
if (j>=0) {
|
|
|
678 |
ierr = QEPCompareEigenvalues(qep,re,0.0,eig[perm[j]],0.0,&result);CHKERRQ(ierr);
|
|
|
679 |
}
|
|
|
680 |
}
|
|
|
681 |
}
|
|
|
682 |
PetscFunctionReturn(0);
|
|
|
683 |
}
|
|
|
684 |
|
|
|
685 |
#undef __FUNCT__
|
|
|
686 |
#define __FUNCT__ "QEPCompareEigenvalues"
|
|
|
687 |
/*@
|
|
|
688 |
QEPCompareEigenvalues - Compares two (possibly complex) eigenvalues according
|
|
|
689 |
to a certain criterion.
|
|
|
690 |
|
|
|
691 |
Not Collective
|
|
|
692 |
|
|
|
693 |
Input Parameters:
|
|
|
694 |
+ qep - the quadratic eigensolver context
|
|
|
695 |
. ar - real part of the 1st eigenvalue
|
|
|
696 |
. ai - imaginary part of the 1st eigenvalue
|
|
|
697 |
. br - real part of the 2nd eigenvalue
|
|
|
698 |
- bi - imaginary part of the 2nd eigenvalue
|
|
|
699 |
|
|
|
700 |
Output Parameter:
|
|
|
701 |
. res - result of comparison
|
|
|
702 |
|
|
|
703 |
Notes:
|
|
|
704 |
Returns an integer less than, equal to, or greater than zero if the first
|
|
|
705 |
eigenvalue is considered to be respectively less than, equal to, or greater
|
|
|
706 |
than the second one.
|
|
|
707 |
|
|
|
708 |
The criterion of comparison is related to the 'which' parameter set with
|
|
|
709 |
QEPSetWhichEigenpairs().
|
|
|
710 |
|
|
|
711 |
Level: developer
|
|
|
712 |
|
|
|
713 |
.seealso: QEPSortEigenvalues(), QEPSetWhichEigenpairs()
|
|
|
714 |
@*/
|
|
|
715 |
PetscErrorCode QEPCompareEigenvalues(QEP qep,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result)
|
|
|
716 |
{
|
|
|
717 |
PetscReal a,b;
|
|
|
718 |
|
|
|
719 |
PetscFunctionBegin;
|
|
|
720 |
switch(qep->which) {
|
|
|
721 |
case QEP_LARGEST_MAGNITUDE:
|
|
|
722 |
case QEP_SMALLEST_MAGNITUDE:
|
|
|
723 |
a = SlepcAbsEigenvalue(ar,ai);
|
|
|
724 |
b = SlepcAbsEigenvalue(br,bi);
|
|
|
725 |
break;
|
|
|
726 |
case QEP_LARGEST_REAL:
|
|
|
727 |
case QEP_SMALLEST_REAL:
|
|
|
728 |
a = PetscRealPart(ar);
|
|
|
729 |
b = PetscRealPart(br);
|
|
|
730 |
break;
|
|
|
731 |
case QEP_LARGEST_IMAGINARY:
|
|
|
732 |
case QEP_SMALLEST_IMAGINARY:
|
|
|
733 |
#if defined(PETSC_USE_COMPLEX)
|
|
|
734 |
a = PetscImaginaryPart(ar);
|
|
|
735 |
b = PetscImaginaryPart(br);
|
|
|
736 |
#else
|
|
|
737 |
a = PetscAbsReal(ai);
|
|
|
738 |
b = PetscAbsReal(bi);
|
|
|
739 |
#endif
|
|
|
740 |
break;
|
|
|
741 |
default: SETERRQ(1,"Wrong value of which");
|
|
|
742 |
}
|
|
|
743 |
switch(qep->which) {
|
|
|
744 |
case QEP_LARGEST_MAGNITUDE:
|
|
|
745 |
case QEP_LARGEST_REAL:
|
|
|
746 |
case QEP_LARGEST_IMAGINARY:
|
|
|
747 |
if (a<b) *result = -1;
|
|
|
748 |
else if (a>b) *result = 1;
|
|
|
749 |
else *result = 0;
|
|
|
750 |
break;
|
|
|
751 |
default:
|
|
|
752 |
if (a>b) *result = -1;
|
|
|
753 |
else if (a<b) *result = 1;
|
|
|
754 |
else *result = 0;
|
|
|
755 |
}
|
|
|
756 |
PetscFunctionReturn(0);
|
|
|
757 |
}
|
|
|
758 |
|
| 1896 |
jroman |
759 |
#undef __FUNCT__
|
|
|
760 |
#define __FUNCT__ "QEPGetOperationCounters"
|
|
|
761 |
/*@
|
|
|
762 |
QEPGetOperationCounters - Gets the total number of matrix-vector products, dot
|
|
|
763 |
products, and linear solve iterations used by the QEP object during the last
|
|
|
764 |
QEPSolve() call.
|
|
|
765 |
|
|
|
766 |
Not Collective
|
|
|
767 |
|
|
|
768 |
Input Parameter:
|
|
|
769 |
. qep - quadratic eigensolver context
|
|
|
770 |
|
|
|
771 |
Output Parameter:
|
|
|
772 |
+ matvecs - number of matrix-vector product operations
|
|
|
773 |
. dots - number of dot product operations
|
|
|
774 |
- lits - number of linear iterations
|
|
|
775 |
|
|
|
776 |
Notes:
|
|
|
777 |
These counters are reset to zero at each successive call to QEPSolve().
|
|
|
778 |
|
|
|
779 |
Level: intermediate
|
|
|
780 |
|
|
|
781 |
@*/
|
|
|
782 |
PetscErrorCode QEPGetOperationCounters(QEP qep,PetscInt* matvecs,PetscInt* dots,PetscInt* lits)
|
|
|
783 |
{
|
|
|
784 |
PetscErrorCode ierr;
|
|
|
785 |
|
|
|
786 |
PetscFunctionBegin;
|
| 2213 |
jroman |
787 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1896 |
jroman |
788 |
if (matvecs) *matvecs = qep->matvecs;
|
|
|
789 |
if (dots) {
|
|
|
790 |
ierr = IPGetOperationCounters(qep->ip,dots);CHKERRQ(ierr);
|
|
|
791 |
}
|
|
|
792 |
if (lits) *lits = qep->linits;
|
|
|
793 |
PetscFunctionReturn(0);
|
|
|
794 |
}
|