Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1887 jroman 1
/*
2
      QEP routines related to the solution process.
3
 
4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2575 eromero 6
   Copyright (c) 2002-2011, Universitat 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
 
2284 jroman 24
#include <private/qepimpl.h>       /*I "slepcqep.h" I*/
1887 jroman 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;
2216 jroman 50
  PetscBool      flg;
1887 jroman 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;
2216 jroman 60
  ierr = PetscOptionsGetBool(((PetscObject)qep)->prefix,"-qep_view_binary",&flg,PETSC_NULL);CHKERRQ(ierr);
1887 jroman 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
 
2381 jroman 70
  if (!qep->ip) { ierr = QEPGetIP(qep,&qep->ip);CHKERRQ(ierr); }
1887 jroman 71
  if (!qep->setupcalled){ ierr = QEPSetUp(qep);CHKERRQ(ierr); }
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;
2313 jroman 75
  ierr = QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,qep->ncv);CHKERRQ(ierr);
1887 jroman 76
 
77
  ierr = PetscLogEventBegin(QEP_Solve,qep,0,0,0);CHKERRQ(ierr);
78
  ierr = (*qep->ops->solve)(qep);CHKERRQ(ierr);
79
  ierr = PetscLogEventEnd(QEP_Solve,qep,0,0,0);CHKERRQ(ierr);
80
 
81
  if (!qep->reason) {
2214 jroman 82
    SETERRQ(((PetscObject)qep)->comm,1,"Internal error, solver returned without setting converged reason");
1887 jroman 83
  }
84
 
2320 jroman 85
#if !defined(PETSC_USE_COMPLEX)
1887 jroman 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
        ierr = VecScale(qep->V[i+1],-1.0);CHKERRQ(ierr);
93
      }
94
      i++;
95
    }
96
  }
97
#endif
98
 
99
  /* sort eigenvalues according to qep->which parameter */
100
  ierr = PetscFree(qep->perm);CHKERRQ(ierr);
101
  if (qep->nconv > 0) {
102
    ierr = PetscMalloc(sizeof(PetscInt)*qep->nconv,&qep->perm);CHKERRQ(ierr);
103
    ierr = QEPSortEigenvalues(qep,qep->nconv,qep->eigr,qep->eigi,qep->perm);CHKERRQ(ierr);
104
  }
105
 
106
  ierr = PetscOptionsGetString(((PetscObject)qep)->prefix,"-qep_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
107
  if (flg && !PetscPreLoadingOn) {
108
    ierr = PetscViewerASCIIOpen(((PetscObject)qep)->comm,filename,&viewer);CHKERRQ(ierr);
109
    ierr = QEPView(qep,viewer);CHKERRQ(ierr);
2305 jroman 110
    ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1887 jroman 111
  }
112
 
113
  flg = PETSC_FALSE;
2216 jroman 114
  ierr = PetscOptionsGetBool(((PetscObject)qep)->prefix,"-qep_plot_eigs",&flg,PETSC_NULL);CHKERRQ(ierr);
1887 jroman 115
  if (flg) {
116
    ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
117
                             PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);CHKERRQ(ierr);
118
    ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
119
    ierr = PetscDrawSPCreate(draw,1,&drawsp);CHKERRQ(ierr);
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
      ierr = PetscDrawSPAddPoint(drawsp,&re,&im);CHKERRQ(ierr);
129
    }
130
    ierr = PetscDrawSPDraw(drawsp);CHKERRQ(ierr);
2305 jroman 131
    ierr = PetscDrawSPDestroy(&drawsp);CHKERRQ(ierr);
132
    ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1887 jroman 133
  }
134
 
2080 eromero 135
  /* Remove the initial subspace */
136
  qep->nini = 0;
1887 jroman 137
  PetscFunctionReturn(0);
138
}
139
 
140
#undef __FUNCT__  
141
#define __FUNCT__ "QEPGetIterationNumber"
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
148
 
149
   Input Parameter:
150
.  qep - the quadratic eigensolver context
151
 
152
   Output Parameter:
153
.  its - number of iterations
154
 
155
   Level: intermediate
156
 
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.
163
 
164
.seealso: QEPGetConvergedReason(), QEPSetTolerances()
165
@*/
166
PetscErrorCode QEPGetIterationNumber(QEP qep,PetscInt *its)
167
{
168
  PetscFunctionBegin;
2213 jroman 169
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 170
  PetscValidIntPointer(its,2);
171
  *its = qep->its;
172
  PetscFunctionReturn(0);
173
}
174
 
175
#undef __FUNCT__  
176
#define __FUNCT__ "QEPGetConverged"
177
/*@
178
   QEPGetConverged - Gets the number of converged eigenpairs.
179
 
180
   Not Collective
181
 
182
   Input Parameter:
183
.  qep - the quadratic eigensolver context
184
 
185
   Output Parameter:
186
.  nconv - number of converged eigenpairs
187
 
188
   Note:
189
   This function should be called after QEPSolve() has finished.
190
 
191
   Level: beginner
192
 
193
.seealso: QEPSetDimensions(), QEPSolve()
194
@*/
195
PetscErrorCode QEPGetConverged(QEP qep,PetscInt *nconv)
196
{
197
  PetscFunctionBegin;
2213 jroman 198
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 199
  PetscValidIntPointer(nconv,2);
200
  *nconv = qep->nconv;
201
  PetscFunctionReturn(0);
202
}
203
 
204
#undef __FUNCT__  
205
#define __FUNCT__ "QEPGetConvergedReason"
206
/*@C
207
   QEPGetConvergedReason - Gets the reason why the QEPSolve() iteration was
208
   stopped.
209
 
210
   Not Collective
211
 
212
   Input Parameter:
213
.  qep - the quadratic eigensolver context
214
 
215
   Output Parameter:
216
.  reason - negative value indicates diverged, positive value converged
217
 
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
222
 
223
   Note:
224
   Can only be called after the call to QEPSolve() is complete.
225
 
226
   Level: intermediate
227
 
228
.seealso: QEPSetTolerances(), QEPSolve(), QEPConvergedReason
229
@*/
230
PetscErrorCode QEPGetConvergedReason(QEP qep,QEPConvergedReason *reason)
231
{
232
  PetscFunctionBegin;
2213 jroman 233
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
2326 jroman 234
  PetscValidPointer(reason,2);
1887 jroman 235
  *reason = qep->reason;
236
  PetscFunctionReturn(0);
237
}
238
 
239
#undef __FUNCT__  
240
#define __FUNCT__ "QEPGetEigenpair"
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.
244
 
245
   Not Collective, but vectors are shared by all processors that share the QEP
246
 
247
   Input Parameters:
248
+  qep - quadratic eigensolver context
249
-  i   - index of the solution
250
 
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
256
 
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).
262
 
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().
266
 
267
   Level: beginner
268
 
269
.seealso: QEPSolve(), QEPGetConverged(), QEPSetWhichEigenpairs()
270
@*/
2326 jroman 271
PetscErrorCode QEPGetEigenpair(QEP qep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
1887 jroman 272
{
273
  PetscInt       k;
274
  PetscErrorCode ierr;
275
 
276
  PetscFunctionBegin;
2213 jroman 277
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
2326 jroman 278
  if (Vr) { PetscValidHeaderSpecific(Vr,VEC_CLASSID,6); PetscCheckSameComm(qep,1,Vr,6); }
279
  if (Vi) { PetscValidHeaderSpecific(Vi,VEC_CLASSID,7); PetscCheckSameComm(qep,1,Vi,7); }
1887 jroman 280
  if (!qep->eigr || !qep->eigi || !qep->V) {
2331 jroman 281
    SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONGSTATE,"QEPSolve must be called first");
1887 jroman 282
  }
283
  if (i<0 || i>=qep->nconv) {
2331 jroman 284
    SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
1887 jroman 285
  }
286
 
287
  if (!qep->perm) k = i;
288
  else k = qep->perm[i];
289
 
290
  /* eigenvalue */
2320 jroman 291
#if defined(PETSC_USE_COMPLEX)
1887 jroman 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 */
2320 jroman 300
#if defined(PETSC_USE_COMPLEX)
1887 jroman 301
  if (Vr) { ierr = VecCopy(qep->V[k],Vr);CHKERRQ(ierr); }
302
  if (Vi) { ierr = VecSet(Vi,0.0);CHKERRQ(ierr); }
303
#else
304
  if (qep->eigi[k]>0) { /* first value of conjugate pair */
305
    if (Vr) { ierr = VecCopy(qep->V[k],Vr);CHKERRQ(ierr); }
306
    if (Vi) { ierr = VecCopy(qep->V[k+1],Vi);CHKERRQ(ierr); }
307
  } else if (qep->eigi[k]<0) { /* second value of conjugate pair */
308
    if (Vr) { ierr = VecCopy(qep->V[k-1],Vr);CHKERRQ(ierr); }
309
    if (Vi) {
310
      ierr = VecCopy(qep->V[k],Vi);CHKERRQ(ierr);
311
      ierr = VecScale(Vi,-1.0);CHKERRQ(ierr);
312
    }
313
  } else { /* real eigenvalue */
314
    if (Vr) { ierr = VecCopy(qep->V[k],Vr);CHKERRQ(ierr); }
315
    if (Vi) { ierr = VecSet(Vi,0.0);CHKERRQ(ierr); }
316
  }
317
#endif
318
  PetscFunctionReturn(0);
319
}
320
 
321
#undef __FUNCT__  
322
#define __FUNCT__ "QEPGetErrorEstimate"
323
/*@
324
   QEPGetErrorEstimate - Returns the error estimate associated to the i-th
325
   computed eigenpair.
326
 
327
   Not Collective
328
 
329
   Input Parameter:
330
+  qep - quadratic eigensolver context
331
-  i   - index of eigenpair
332
 
333
   Output Parameter:
334
.  errest - the error estimate
335
 
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.
340
 
341
   Level: advanced
342
 
343
.seealso: QEPComputeRelativeError()
344
@*/
2326 jroman 345
PetscErrorCode QEPGetErrorEstimate(QEP qep,PetscInt i,PetscReal *errest)
1887 jroman 346
{
347
  PetscFunctionBegin;
2213 jroman 348
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
2326 jroman 349
  PetscValidPointer(errest,3);
1887 jroman 350
  if (!qep->eigr || !qep->eigi) {
2331 jroman 351
    SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONGSTATE,"QEPSolve must be called first");
1887 jroman 352
  }
353
  if (i<0 || i>=qep->nconv) {
2331 jroman 354
    SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
1887 jroman 355
  }
356
  if (qep->perm) i = qep->perm[i];  
357
  if (errest) *errest = qep->errest[i];
358
  PetscFunctionReturn(0);
359
}
360
 
361
#undef __FUNCT__  
362
#define __FUNCT__ "QEPComputeResidualNorm_Private"
363
/*
364
   QEPComputeResidualNorm_Private - Computes the norm of the residual vector
365
   associated with an eigenpair.
366
*/
2331 jroman 367
PetscErrorCode QEPComputeResidualNorm_Private(QEP qep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *norm)
1887 jroman 368
{
369
  PetscErrorCode ierr;
370
  Vec            u,w;
1888 jroman 371
  Mat            M=qep->M,C=qep->C,K=qep->K;
2320 jroman 372
#if !defined(PETSC_USE_COMPLEX)
1887 jroman 373
  Vec            v,y,z;
374
  PetscReal      ni,nr;
375
  PetscScalar    a1,a2;
376
#endif
377
 
378
  PetscFunctionBegin;
1952 jroman 379
  ierr = VecDuplicate(qep->V[0],&u);CHKERRQ(ierr);
1887 jroman 380
  ierr = VecDuplicate(u,&w);CHKERRQ(ierr);
381
 
2320 jroman 382
#if !defined(PETSC_USE_COMPLEX)
1887 jroman 383
  if (ki == 0 ||
384
    PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
385
#endif
1888 jroman 386
    ierr = MatMult(K,xr,u);CHKERRQ(ierr);                 /* u=K*x */
1887 jroman 387
    if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
388
      ierr = MatMult(C,xr,w);CHKERRQ(ierr);               /* w=C*x */
1888 jroman 389
      ierr = VecAXPY(u,kr,w);CHKERRQ(ierr);               /* u=l*C*x+K*x */
1893 jroman 390
      ierr = MatMult(M,xr,w);CHKERRQ(ierr);               /* w=M*x */
1888 jroman 391
      ierr = VecAXPY(u,kr*kr,w);CHKERRQ(ierr);            /* u=l^2*M*x+l*C*x+K*x */
1887 jroman 392
    }
393
    ierr = VecNorm(u,NORM_2,norm);CHKERRQ(ierr);  
2320 jroman 394
#if !defined(PETSC_USE_COMPLEX)
1887 jroman 395
  } else {
396
    ierr = VecDuplicate(u,&v);CHKERRQ(ierr);
397
    ierr = VecDuplicate(u,&y);CHKERRQ(ierr);
398
    ierr = VecDuplicate(u,&z);CHKERRQ(ierr);
399
    a1 = kr*kr-ki*ki;
400
    a2 = 2.0*kr*ki;
1893 jroman 401
    ierr = MatMult(K,xr,u);CHKERRQ(ierr);           /* u=K*xr */
402
    if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
403
      ierr = MatMult(C,xr,v);CHKERRQ(ierr);         /* v=C*xr */
404
      ierr = MatMult(C,xi,w);CHKERRQ(ierr);         /* w=C*xi */
405
      ierr = MatMult(M,xr,y);CHKERRQ(ierr);         /* y=M*xr */
406
      ierr = MatMult(M,xi,z);CHKERRQ(ierr);         /* z=M*xi */
407
      ierr = VecAXPY(u,kr,v);CHKERRQ(ierr);         /* u=kr*C*xr+K*xr */
408
      ierr = VecAXPY(u,-ki,w);CHKERRQ(ierr);        /* u=kr*C*xr-ki*C*xi+K*xr */
409
      ierr = VecAXPY(u,a1,y);CHKERRQ(ierr);         /* u=a1*M*xr+kr*C*xr-ki*C*xi+K*xr */
410
      ierr = VecAXPY(u,-a2,z);CHKERRQ(ierr);        /* u=a1*M*xr-a2*M*xi+kr*C*xr-ki*C*xi+K*xr */
411
    }
1887 jroman 412
    ierr = VecNorm(u,NORM_2,&nr);CHKERRQ(ierr);
1888 jroman 413
    ierr = MatMult(K,xi,u);CHKERRQ(ierr);         /* u=K*xi */
1893 jroman 414
    if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
415
      ierr = VecAXPY(u,kr,w);CHKERRQ(ierr);         /* u=kr*C*xi+K*xi */
416
      ierr = VecAXPY(u,ki,v);CHKERRQ(ierr);         /* u=kr*C*xi+ki*C*xi+K*xi */
417
      ierr = VecAXPY(u,a1,z);CHKERRQ(ierr);         /* u=a1*M*xi+kr*C*xi+ki*C*xi+K*xi */
418
      ierr = VecAXPY(u,a2,y);CHKERRQ(ierr);         /* u=a1*M*xi+a2*M*ki+kr*C*xi+ki*C*xi+K*xi */
419
    }
1887 jroman 420
    ierr = VecNorm(u,NORM_2,&ni);CHKERRQ(ierr);
421
    *norm = SlepcAbsEigenvalue(nr,ni);
2305 jroman 422
    ierr = VecDestroy(&v);CHKERRQ(ierr);
423
    ierr = VecDestroy(&y);CHKERRQ(ierr);
424
    ierr = VecDestroy(&z);CHKERRQ(ierr);
1887 jroman 425
  }
426
#endif
427
 
2305 jroman 428
  ierr = VecDestroy(&w);CHKERRQ(ierr);
429
  ierr = VecDestroy(&u);CHKERRQ(ierr);
1887 jroman 430
  PetscFunctionReturn(0);
431
}
432
 
433
#undef __FUNCT__  
434
#define __FUNCT__ "QEPComputeResidualNorm"
435
/*@
436
   QEPComputeResidualNorm - Computes the norm of the residual vector associated with
437
   the i-th computed eigenpair.
438
 
439
   Collective on QEP
440
 
441
   Input Parameter:
442
.  qep - the quadratic eigensolver context
443
.  i   - the solution index
444
 
445
   Output Parameter:
1888 jroman 446
.  norm - the residual norm, computed as ||(l^2*M+l*C+K)x||_2 where l is the
1887 jroman 447
   eigenvalue and x is the eigenvector.
1893 jroman 448
   If l=0 then the residual norm is computed as ||Kx||_2.
1887 jroman 449
 
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().
454
 
455
   Level: beginner
456
 
457
.seealso: QEPSolve(), QEPGetConverged(), QEPSetWhichEigenpairs()
458
@*/
2331 jroman 459
PetscErrorCode QEPComputeResidualNorm(QEP qep,PetscInt i,PetscReal *norm)
1887 jroman 460
{
461
  PetscErrorCode ierr;
462
  Vec            xr,xi;
463
  PetscScalar    kr,ki;
464
 
465
  PetscFunctionBegin;
2213 jroman 466
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
2328 jroman 467
  PetscValidLogicalCollectiveInt(qep,i,2);
1887 jroman 468
  PetscValidPointer(norm,3);
1952 jroman 469
  ierr = VecDuplicate(qep->V[0],&xr);CHKERRQ(ierr);
470
  ierr = VecDuplicate(qep->V[0],&xi);CHKERRQ(ierr);
1887 jroman 471
  ierr = QEPGetEigenpair(qep,i,&kr,&ki,xr,xi);CHKERRQ(ierr);
472
  ierr = QEPComputeResidualNorm_Private(qep,kr,ki,xr,xi,norm);CHKERRQ(ierr);
2305 jroman 473
  ierr = VecDestroy(&xr);CHKERRQ(ierr);
474
  ierr = VecDestroy(&xi);CHKERRQ(ierr);
1887 jroman 475
  PetscFunctionReturn(0);
476
}
477
 
478
#undef __FUNCT__  
479
#define __FUNCT__ "QEPComputeRelativeError_Private"
480
/*
481
   QEPComputeRelativeError_Private - Computes the relative error bound
482
   associated with an eigenpair.
483
*/
2331 jroman 484
PetscErrorCode QEPComputeRelativeError_Private(QEP qep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *error)
1887 jroman 485
{
486
  PetscErrorCode ierr;
2331 jroman 487
  PetscReal      norm,er;
2320 jroman 488
#if !defined(PETSC_USE_COMPLEX)
1887 jroman 489
  PetscReal      ei;
490
#endif
491
 
492
  PetscFunctionBegin;
1893 jroman 493
  ierr = QEPComputeResidualNorm_Private(qep,kr,ki,xr,xi,&norm);CHKERRQ(ierr);
2320 jroman 494
#if !defined(PETSC_USE_COMPLEX)
1887 jroman 495
  if (ki == 0 ||
496
    PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
497
#endif
1893 jroman 498
    ierr = VecNorm(xr,NORM_2,&er);CHKERRQ(ierr);
1887 jroman 499
    if (PetscAbsScalar(kr) > norm) {
1893 jroman 500
      *error = norm/(PetscAbsScalar(kr)*er);
1887 jroman 501
    } else {
1893 jroman 502
      *error = norm/er;
1887 jroman 503
    }
2320 jroman 504
#if !defined(PETSC_USE_COMPLEX)
1887 jroman 505
  } else {
1893 jroman 506
    ierr = VecNorm(xr,NORM_2,&er);CHKERRQ(ierr);  
507
    ierr = VecNorm(xi,NORM_2,&ei);CHKERRQ(ierr);  
1887 jroman 508
    if (SlepcAbsEigenvalue(kr,ki) > norm) {
1893 jroman 509
      *error = norm/(SlepcAbsEigenvalue(kr,ki)*SlepcAbsEigenvalue(er,ei));
1887 jroman 510
    } else {
1893 jroman 511
      *error = norm/SlepcAbsEigenvalue(er,ei);
1887 jroman 512
    }
513
  }
514
#endif    
515
  PetscFunctionReturn(0);
516
}
517
 
518
#undef __FUNCT__  
519
#define __FUNCT__ "QEPComputeRelativeError"
520
/*@
521
   QEPComputeRelativeError - Computes the relative error bound associated
522
   with the i-th computed eigenpair.
523
 
524
   Collective on QEP
525
 
526
   Input Parameter:
527
.  qep - the quadratic eigensolver context
528
.  i   - the solution index
529
 
530
   Output Parameter:
1888 jroman 531
.  error - the relative error bound, computed as ||(l^2*M+l*C+K)x||_2/||lx||_2 where
1887 jroman 532
   l is the eigenvalue and x is the eigenvector.
1893 jroman 533
   If l=0 the relative error is computed as ||Kx||_2/||x||_2.
1887 jroman 534
 
535
   Level: beginner
536
 
537
.seealso: QEPSolve(), QEPComputeResidualNorm(), QEPGetErrorEstimate()
538
@*/
2331 jroman 539
PetscErrorCode QEPComputeRelativeError(QEP qep,PetscInt i,PetscReal *error)
1887 jroman 540
{
541
  PetscErrorCode ierr;
542
  Vec            xr,xi;  
543
  PetscScalar    kr,ki;  
544
 
545
  PetscFunctionBegin;
2213 jroman 546
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);  
2328 jroman 547
  PetscValidLogicalCollectiveInt(qep,i,2);
1887 jroman 548
  PetscValidPointer(error,3);
1952 jroman 549
  ierr = VecDuplicate(qep->V[0],&xr);CHKERRQ(ierr);
550
  ierr = VecDuplicate(qep->V[0],&xi);CHKERRQ(ierr);
1887 jroman 551
  ierr = QEPGetEigenpair(qep,i,&kr,&ki,xr,xi);CHKERRQ(ierr);
552
  ierr = QEPComputeRelativeError_Private(qep,kr,ki,xr,xi,error);CHKERRQ(ierr);  
2305 jroman 553
  ierr = VecDestroy(&xr);CHKERRQ(ierr);
554
  ierr = VecDestroy(&xi);CHKERRQ(ierr);
1887 jroman 555
  PetscFunctionReturn(0);
556
}
557
 
558
#undef __FUNCT__  
559
#define __FUNCT__ "QEPSortEigenvalues"
560
/*@
561
   QEPSortEigenvalues - Sorts a list of eigenvalues according to the criterion
562
   specified via QEPSetWhichEigenpairs().
563
 
564
   Not Collective
565
 
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)
571
 
572
   Output Parameter:
573
.  perm  - resulting permutation
574
 
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.
579
 
580
   Level: developer
581
 
582
.seealso: QEPSortEigenvaluesReal(), QEPSetWhichEigenpairs()
583
@*/
584
PetscErrorCode QEPSortEigenvalues(QEP qep,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm)
585
{
586
  PetscErrorCode ierr;
587
  PetscScalar    re,im;
588
  PetscInt       i,j,result,tmp;
589
 
590
  PetscFunctionBegin;
2349 jroman 591
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);  
592
  PetscValidScalarPointer(eigr,3);
593
  PetscValidScalarPointer(eigi,4);
594
  PetscValidIntPointer(perm,5);
1887 jroman 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;
2320 jroman 601
#if !defined(PETSC_USE_COMPLEX)
1887 jroman 602
    if (im != 0) {
603
      /* complex eigenvalue */
604
      i--;
605
      im = eigi[perm[i]];
606
    }
607
#endif
608
    while (j<n) {
609
      ierr = QEPCompareEigenvalues(qep,re,im,eigr[perm[j]],eigi[perm[j]],&result);CHKERRQ(ierr);
610
      if (result >= 0) break;
2320 jroman 611
#if !defined(PETSC_USE_COMPLEX)
1887 jroman 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++;
2320 jroman 618
#if !defined(PETSC_USE_COMPLEX)
1887 jroman 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
  PetscFunctionReturn(0);
637
}
638
 
639
#undef __FUNCT__  
640
#define __FUNCT__ "QEPSortEigenvaluesReal"
641
/*@
642
   QEPSortEigenvaluesReal - Sorts a list of eigenvalues according to a certain
643
   criterion (version for real eigenvalues only).
644
 
645
   Not Collective
646
 
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)
651
 
652
   Output Parameter:
653
.  perm  - resulting permutation
654
 
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.
659
 
660
   Level: developer
661
 
662
.seealso: QEPSortEigenvalues(), QEPSetWhichEigenpairs(), QEPCompareEigenvalues()
663
@*/
664
PetscErrorCode QEPSortEigenvaluesReal(QEP qep,PetscInt n,PetscReal *eig,PetscInt *perm)
665
{
666
  PetscErrorCode ierr;
667
  PetscScalar    re;
668
  PetscInt       i,j,result,tmp;
669
 
670
  PetscFunctionBegin;
2349 jroman 671
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);  
672
  PetscValidScalarPointer(eig,3);
673
  PetscValidIntPointer(perm,4);
1887 jroman 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
    ierr = QEPCompareEigenvalues(qep,re,0.0,eig[perm[j]],0.0,&result);CHKERRQ(ierr);
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
        ierr = QEPCompareEigenvalues(qep,re,0.0,eig[perm[j]],0.0,&result);CHKERRQ(ierr);
684
      }
685
    }
686
  }
687
  PetscFunctionReturn(0);
688
}
689
 
690
#undef __FUNCT__  
691
#define __FUNCT__ "QEPCompareEigenvalues"
692
/*@
693
   QEPCompareEigenvalues - Compares two (possibly complex) eigenvalues according
694
   to a certain criterion.
695
 
696
   Not Collective
697
 
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
704
 
705
   Output Parameter:
706
.  res    - result of comparison
707
 
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.
712
 
713
   The criterion of comparison is related to the 'which' parameter set with
714
   QEPSetWhichEigenpairs().
715
 
716
   Level: developer
717
 
718
.seealso: QEPSortEigenvalues(), QEPSetWhichEigenpairs()
719
@*/
720
PetscErrorCode QEPCompareEigenvalues(QEP qep,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result)
721
{
2317 jroman 722
  PetscReal a,b;
1887 jroman 723
 
724
  PetscFunctionBegin;
2349 jroman 725
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);  
726
  PetscValidIntPointer(result,6);
1887 jroman 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;
2214 jroman 748
    default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of which");
1887 jroman 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
  PetscFunctionReturn(0);
764
}
765
 
1896 jroman 766
#undef __FUNCT__  
767
#define __FUNCT__ "QEPGetOperationCounters"
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.
772
 
773
   Not Collective
774
 
775
   Input Parameter:
776
.  qep - quadratic eigensolver context
777
 
778
   Output Parameter:
779
+  matvecs - number of matrix-vector product operations
780
.  dots    - number of dot product operations
781
-  lits    - number of linear iterations
782
 
783
   Notes:
784
   These counters are reset to zero at each successive call to QEPSolve().
785
 
786
   Level: intermediate
787
 
788
@*/
789
PetscErrorCode QEPGetOperationCounters(QEP qep,PetscInt* matvecs,PetscInt* dots,PetscInt* lits)
790
{
791
  PetscErrorCode ierr;
792
 
793
  PetscFunctionBegin;
2213 jroman 794
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1896 jroman 795
  if (matvecs) *matvecs = qep->matvecs;
796
  if (dots) {
2381 jroman 797
    if (!qep->ip) { ierr = QEPGetIP(qep,&qep->ip);CHKERRQ(ierr); }
1896 jroman 798
    ierr = IPGetOperationCounters(qep->ip,dots);CHKERRQ(ierr);
799
  }
800
  if (lits) *lits = qep->linits;
801
  PetscFunctionReturn(0);
802
}