Subversion Repositories slepc-dev

Rev

Go to most recent revision | 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
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
 
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
 
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;
2313 jroman 77
  ierr = QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,qep->ncv);CHKERRQ(ierr);
1887 jroman 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) {
2214 jroman 84
    SETERRQ(((PetscObject)qep)->comm,1,"Internal error, solver returned without setting converged reason");
1887 jroman 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);
2305 jroman 112
    ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1887 jroman 113
  }
114
 
115
  flg = PETSC_FALSE;
2216 jroman 116
  ierr = PetscOptionsGetBool(((PetscObject)qep)->prefix,"-qep_plot_eigs",&flg,PETSC_NULL);CHKERRQ(ierr);
1887 jroman 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);
2305 jroman 133
    ierr = PetscDrawSPDestroy(&drawsp);CHKERRQ(ierr);
134
    ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1887 jroman 135
  }
136
 
2080 eromero 137
  /* Remove the initial subspace */
138
  qep->nini = 0;
1887 jroman 139
  PetscFunctionReturn(0);
140
}
141
 
142
#undef __FUNCT__  
143
#define __FUNCT__ "QEPGetIterationNumber"
144
/*@
145
   QEPGetIterationNumber - Gets the current iteration number. If the
146
   call to QEPSolve() is complete, then it returns the number of iterations
147
   carried out by the solution method.
148
 
149
   Not Collective
150
 
151
   Input Parameter:
152
.  qep - the quadratic eigensolver context
153
 
154
   Output Parameter:
155
.  its - number of iterations
156
 
157
   Level: intermediate
158
 
159
   Note:
160
   During the i-th iteration this call returns i-1. If QEPSolve() is
161
   complete, then parameter "its" contains either the iteration number at
162
   which convergence was successfully reached, or failure was detected.  
163
   Call QEPGetConvergedReason() to determine if the solver converged or
164
   failed and why.
165
 
166
.seealso: QEPGetConvergedReason(), QEPSetTolerances()
167
@*/
168
PetscErrorCode QEPGetIterationNumber(QEP qep,PetscInt *its)
169
{
170
  PetscFunctionBegin;
2213 jroman 171
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 172
  PetscValidIntPointer(its,2);
173
  *its = qep->its;
174
  PetscFunctionReturn(0);
175
}
176
 
177
#undef __FUNCT__  
178
#define __FUNCT__ "QEPGetConverged"
179
/*@
180
   QEPGetConverged - Gets the number of converged eigenpairs.
181
 
182
   Not Collective
183
 
184
   Input Parameter:
185
.  qep - the quadratic eigensolver context
186
 
187
   Output Parameter:
188
.  nconv - number of converged eigenpairs
189
 
190
   Note:
191
   This function should be called after QEPSolve() has finished.
192
 
193
   Level: beginner
194
 
195
.seealso: QEPSetDimensions(), QEPSolve()
196
@*/
197
PetscErrorCode QEPGetConverged(QEP qep,PetscInt *nconv)
198
{
199
  PetscFunctionBegin;
2213 jroman 200
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 201
  PetscValidIntPointer(nconv,2);
202
  *nconv = qep->nconv;
203
  PetscFunctionReturn(0);
204
}
205
 
206
 
207
#undef __FUNCT__  
208
#define __FUNCT__ "QEPGetConvergedReason"
209
/*@C
210
   QEPGetConvergedReason - Gets the reason why the QEPSolve() iteration was
211
   stopped.
212
 
213
   Not Collective
214
 
215
   Input Parameter:
216
.  qep - the quadratic eigensolver context
217
 
218
   Output Parameter:
219
.  reason - negative value indicates diverged, positive value converged
220
 
221
   Possible values for reason:
222
+  QEP_CONVERGED_TOL - converged up to tolerance
223
.  QEP_DIVERGED_ITS - required more than its to reach convergence
224
-  QEP_DIVERGED_BREAKDOWN - generic breakdown in method
225
 
226
   Note:
227
   Can only be called after the call to QEPSolve() is complete.
228
 
229
   Level: intermediate
230
 
231
.seealso: QEPSetTolerances(), QEPSolve(), QEPConvergedReason
232
@*/
233
PetscErrorCode QEPGetConvergedReason(QEP qep,QEPConvergedReason *reason)
234
{
235
  PetscFunctionBegin;
2213 jroman 236
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 237
  PetscValidIntPointer(reason,2);
238
  *reason = qep->reason;
239
  PetscFunctionReturn(0);
240
}
241
 
242
#undef __FUNCT__  
243
#define __FUNCT__ "QEPGetEigenpair"
244
/*@
245
   QEPGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
246
   QEPSolve(). The solution consists in both the eigenvalue and the eigenvector.
247
 
248
   Not Collective, but vectors are shared by all processors that share the QEP
249
 
250
   Input Parameters:
251
+  qep - quadratic eigensolver context
252
-  i   - index of the solution
253
 
254
   Output Parameters:
255
+  eigr - real part of eigenvalue
256
.  eigi - imaginary part of eigenvalue
257
.  Vr   - real part of eigenvector
258
-  Vi   - imaginary part of eigenvector
259
 
260
   Notes:
261
   If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
262
   configured with complex scalars the eigenvalue is stored
263
   directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
264
   set to zero).
265
 
266
   The index i should be a value between 0 and nconv-1 (see QEPGetConverged()).
267
   Eigenpairs are indexed according to the ordering criterion established
268
   with QEPSetWhichEigenpairs().
269
 
270
   Level: beginner
271
 
272
.seealso: QEPSolve(), QEPGetConverged(), QEPSetWhichEigenpairs()
273
@*/
274
PetscErrorCode QEPGetEigenpair(QEP qep, PetscInt i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
275
{
276
  PetscInt       k;
277
  PetscErrorCode ierr;
278
 
279
  PetscFunctionBegin;
2213 jroman 280
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 281
  if (!qep->eigr || !qep->eigi || !qep->V) {
2214 jroman 282
    SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONGSTATE, "QEPSolve must be called first");
1887 jroman 283
  }
284
  if (i<0 || i>=qep->nconv) {
2214 jroman 285
    SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
1887 jroman 286
  }
287
 
288
  if (!qep->perm) k = i;
289
  else k = qep->perm[i];
290
 
291
  /* eigenvalue */
292
#ifdef PETSC_USE_COMPLEX
293
  if (eigr) *eigr = qep->eigr[k];
294
  if (eigi) *eigi = 0;
295
#else
296
  if (eigr) *eigr = qep->eigr[k];
297
  if (eigi) *eigi = qep->eigi[k];
298
#endif
299
 
300
  /* eigenvector */
301
#ifdef PETSC_USE_COMPLEX
302
  if (Vr) { ierr = VecCopy(qep->V[k],Vr);CHKERRQ(ierr); }
303
  if (Vi) { ierr = VecSet(Vi,0.0);CHKERRQ(ierr); }
304
#else
305
  if (qep->eigi[k]>0) { /* first value of conjugate pair */
306
    if (Vr) { ierr = VecCopy(qep->V[k],Vr);CHKERRQ(ierr); }
307
    if (Vi) { ierr = VecCopy(qep->V[k+1],Vi);CHKERRQ(ierr); }
308
  } else if (qep->eigi[k]<0) { /* second value of conjugate pair */
309
    if (Vr) { ierr = VecCopy(qep->V[k-1],Vr);CHKERRQ(ierr); }
310
    if (Vi) {
311
      ierr = VecCopy(qep->V[k],Vi);CHKERRQ(ierr);
312
      ierr = VecScale(Vi,-1.0);CHKERRQ(ierr);
313
    }
314
  } else { /* real eigenvalue */
315
    if (Vr) { ierr = VecCopy(qep->V[k],Vr);CHKERRQ(ierr); }
316
    if (Vi) { ierr = VecSet(Vi,0.0);CHKERRQ(ierr); }
317
  }
318
#endif
319
  PetscFunctionReturn(0);
320
}
321
 
322
#undef __FUNCT__  
323
#define __FUNCT__ "QEPGetErrorEstimate"
324
/*@
325
   QEPGetErrorEstimate - Returns the error estimate associated to the i-th
326
   computed eigenpair.
327
 
328
   Not Collective
329
 
330
   Input Parameter:
331
+  qep - quadratic eigensolver context
332
-  i   - index of eigenpair
333
 
334
   Output Parameter:
335
.  errest - the error estimate
336
 
337
   Notes:
338
   This is the error estimate used internally by the eigensolver. The actual
339
   error bound can be computed with QEPComputeRelativeError(). See also the users
340
   manual for details.
341
 
342
   Level: advanced
343
 
344
.seealso: QEPComputeRelativeError()
345
@*/
346
PetscErrorCode QEPGetErrorEstimate(QEP qep, PetscInt i, PetscReal *errest)
347
{
348
  PetscFunctionBegin;
2213 jroman 349
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 350
  if (!qep->eigr || !qep->eigi) {
2214 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) {
2214 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
*/
367
PetscErrorCode QEPComputeResidualNorm_Private(QEP qep, PetscScalar kr, PetscScalar ki, Vec xr, Vec xi, PetscReal *norm)
368
{
369
  PetscErrorCode ierr;
370
  Vec            u,w;
1888 jroman 371
  Mat            M=qep->M,C=qep->C,K=qep->K;
1887 jroman 372
#ifndef PETSC_USE_COMPLEX
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
 
382
#ifndef PETSC_USE_COMPLEX
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);  
394
#ifndef PETSC_USE_COMPLEX
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
@*/
459
PetscErrorCode QEPComputeResidualNorm(QEP qep, PetscInt i, PetscReal *norm)
460
{
461
  PetscErrorCode ierr;
462
  Vec            xr,xi;
463
  PetscScalar    kr,ki;
464
 
465
  PetscFunctionBegin;
2213 jroman 466
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 467
  PetscValidPointer(norm,3);
1952 jroman 468
  ierr = VecDuplicate(qep->V[0],&xr);CHKERRQ(ierr);
469
  ierr = VecDuplicate(qep->V[0],&xi);CHKERRQ(ierr);
1887 jroman 470
  ierr = QEPGetEigenpair(qep,i,&kr,&ki,xr,xi);CHKERRQ(ierr);
471
  ierr = QEPComputeResidualNorm_Private(qep,kr,ki,xr,xi,norm);CHKERRQ(ierr);
2305 jroman 472
  ierr = VecDestroy(&xr);CHKERRQ(ierr);
473
  ierr = VecDestroy(&xi);CHKERRQ(ierr);
1887 jroman 474
  PetscFunctionReturn(0);
475
}
476
 
477
#undef __FUNCT__  
478
#define __FUNCT__ "QEPComputeRelativeError_Private"
479
/*
480
   QEPComputeRelativeError_Private - Computes the relative error bound
481
   associated with an eigenpair.
482
*/
483
PetscErrorCode QEPComputeRelativeError_Private(QEP qep, PetscScalar kr, PetscScalar ki, Vec xr, Vec xi, PetscReal *error)
484
{
485
  PetscErrorCode ierr;
486
  PetscReal      norm, er;
487
#ifndef PETSC_USE_COMPLEX
488
  PetscReal      ei;
489
#endif
490
 
491
  PetscFunctionBegin;
1893 jroman 492
  ierr = QEPComputeResidualNorm_Private(qep,kr,ki,xr,xi,&norm);CHKERRQ(ierr);
1887 jroman 493
#ifndef PETSC_USE_COMPLEX
494
  if (ki == 0 ||
495
    PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
496
#endif
1893 jroman 497
    ierr = VecNorm(xr,NORM_2,&er);CHKERRQ(ierr);
1887 jroman 498
    if (PetscAbsScalar(kr) > norm) {
1893 jroman 499
      *error = norm/(PetscAbsScalar(kr)*er);
1887 jroman 500
    } else {
1893 jroman 501
      *error = norm/er;
1887 jroman 502
    }
503
#ifndef PETSC_USE_COMPLEX
504
  } else {
1893 jroman 505
    ierr = VecNorm(xr,NORM_2,&er);CHKERRQ(ierr);  
506
    ierr = VecNorm(xi,NORM_2,&ei);CHKERRQ(ierr);  
1887 jroman 507
    if (SlepcAbsEigenvalue(kr,ki) > norm) {
1893 jroman 508
      *error = norm/(SlepcAbsEigenvalue(kr,ki)*SlepcAbsEigenvalue(er,ei));
1887 jroman 509
    } else {
1893 jroman 510
      *error = norm/SlepcAbsEigenvalue(er,ei);
1887 jroman 511
    }
512
  }
513
#endif    
514
  PetscFunctionReturn(0);
515
}
516
 
517
#undef __FUNCT__  
518
#define __FUNCT__ "QEPComputeRelativeError"
519
/*@
520
   QEPComputeRelativeError - Computes the relative error bound associated
521
   with the i-th computed eigenpair.
522
 
523
   Collective on QEP
524
 
525
   Input Parameter:
526
.  qep - the quadratic eigensolver context
527
.  i   - the solution index
528
 
529
   Output Parameter:
1888 jroman 530
.  error - the relative error bound, computed as ||(l^2*M+l*C+K)x||_2/||lx||_2 where
1887 jroman 531
   l is the eigenvalue and x is the eigenvector.
1893 jroman 532
   If l=0 the relative error is computed as ||Kx||_2/||x||_2.
1887 jroman 533
 
534
   Level: beginner
535
 
536
.seealso: QEPSolve(), QEPComputeResidualNorm(), QEPGetErrorEstimate()
537
@*/
538
PetscErrorCode QEPComputeRelativeError(QEP qep, PetscInt i, PetscReal *error)
539
{
540
  PetscErrorCode ierr;
541
  Vec            xr,xi;  
542
  PetscScalar    kr,ki;  
543
 
544
  PetscFunctionBegin;
2213 jroman 545
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);  
1887 jroman 546
  PetscValidPointer(error,3);
1952 jroman 547
  ierr = VecDuplicate(qep->V[0],&xr);CHKERRQ(ierr);
548
  ierr = VecDuplicate(qep->V[0],&xi);CHKERRQ(ierr);
1887 jroman 549
  ierr = QEPGetEigenpair(qep,i,&kr,&ki,xr,xi);CHKERRQ(ierr);
550
  ierr = QEPComputeRelativeError_Private(qep,kr,ki,xr,xi,error);CHKERRQ(ierr);  
2305 jroman 551
  ierr = VecDestroy(&xr);CHKERRQ(ierr);
552
  ierr = VecDestroy(&xi);CHKERRQ(ierr);
1887 jroman 553
  PetscFunctionReturn(0);
554
}
555
 
556
#undef __FUNCT__  
557
#define __FUNCT__ "QEPSortEigenvalues"
558
/*@
559
   QEPSortEigenvalues - Sorts a list of eigenvalues according to the criterion
560
   specified via QEPSetWhichEigenpairs().
561
 
562
   Not Collective
563
 
564
   Input Parameters:
565
+  qep   - the quadratic eigensolver context
566
.  n     - number of eigenvalues in the list
567
.  eigr  - pointer to the array containing the eigenvalues
568
-  eigi  - imaginary part of the eigenvalues (only when using real numbers)
569
 
570
   Output Parameter:
571
.  perm  - resulting permutation
572
 
573
   Note:
574
   The result is a list of indices in the original eigenvalue array
575
   corresponding to the first nev eigenvalues sorted in the specified
576
   criterion.
577
 
578
   Level: developer
579
 
580
.seealso: QEPSortEigenvaluesReal(), QEPSetWhichEigenpairs()
581
@*/
582
PetscErrorCode QEPSortEigenvalues(QEP qep,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm)
583
{
584
  PetscErrorCode ierr;
585
  PetscScalar    re,im;
586
  PetscInt       i,j,result,tmp;
587
 
588
  PetscFunctionBegin;
589
  for (i=0; i<n; i++) { perm[i] = i; }
590
  /* insertion sort */
591
  for (i=n-1; i>=0; i--) {
592
    re = eigr[perm[i]];
593
    im = eigi[perm[i]];
594
    j = i + 1;
595
#ifndef PETSC_USE_COMPLEX
596
    if (im != 0) {
597
      /* complex eigenvalue */
598
      i--;
599
      im = eigi[perm[i]];
600
    }
601
#endif
602
    while (j<n) {
603
      ierr = QEPCompareEigenvalues(qep,re,im,eigr[perm[j]],eigi[perm[j]],&result);CHKERRQ(ierr);
604
      if (result >= 0) break;
605
#ifndef PETSC_USE_COMPLEX
606
      /* keep together every complex conjugated eigenpair */
607
      if (im == 0) {
608
        if (eigi[perm[j]] == 0) {
609
#endif
610
          tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = tmp;
611
          j++;
612
#ifndef PETSC_USE_COMPLEX
613
        } else {
614
          tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp;
615
          j+=2;
616
        }
617
      } else {
618
        if (eigi[perm[j]] == 0) {
619
          tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = perm[j-1]; perm[j-1] = tmp;
620
          j++;
621
        } else {
622
          tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = tmp;
623
          tmp = perm[j-1]; perm[j-1] = perm[j+1]; perm[j+1] = tmp;
624
          j+=2;
625
        }
626
      }
627
#endif
628
    }
629
  }
630
  PetscFunctionReturn(0);
631
}
632
 
633
#undef __FUNCT__  
634
#define __FUNCT__ "QEPSortEigenvaluesReal"
635
/*@
636
   QEPSortEigenvaluesReal - Sorts a list of eigenvalues according to a certain
637
   criterion (version for real eigenvalues only).
638
 
639
   Not Collective
640
 
641
   Input Parameters:
642
+  qep   - the quadratic eigensolver context
643
.  n     - number of eigenvalue in the list
644
-  eig   - pointer to the array containing the eigenvalues (real)
645
 
646
   Output Parameter:
647
.  perm  - resulting permutation
648
 
649
   Note:
650
   The result is a list of indices in the original eigenvalue array
651
   corresponding to the first nev eigenvalues sorted in the specified
652
   criterion.
653
 
654
   Level: developer
655
 
656
.seealso: QEPSortEigenvalues(), QEPSetWhichEigenpairs(), QEPCompareEigenvalues()
657
@*/
658
PetscErrorCode QEPSortEigenvaluesReal(QEP qep,PetscInt n,PetscReal *eig,PetscInt *perm)
659
{
660
  PetscErrorCode ierr;
661
  PetscScalar    re;
662
  PetscInt       i,j,result,tmp;
663
 
664
  PetscFunctionBegin;
665
  for (i=0; i<n; i++) { perm[i] = i; }
666
  /* insertion sort */
667
  for (i=1; i<n; i++) {
668
    re = eig[perm[i]];
669
    j = i-1;
670
    ierr = QEPCompareEigenvalues(qep,re,0.0,eig[perm[j]],0.0,&result);CHKERRQ(ierr);
671
    while (result>0 && j>=0) {
672
      tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
673
      if (j>=0) {
674
        ierr = QEPCompareEigenvalues(qep,re,0.0,eig[perm[j]],0.0,&result);CHKERRQ(ierr);
675
      }
676
    }
677
  }
678
  PetscFunctionReturn(0);
679
}
680
 
681
#undef __FUNCT__  
682
#define __FUNCT__ "QEPCompareEigenvalues"
683
/*@
684
   QEPCompareEigenvalues - Compares two (possibly complex) eigenvalues according
685
   to a certain criterion.
686
 
687
   Not Collective
688
 
689
   Input Parameters:
690
+  qep    - the quadratic eigensolver context
691
.  ar     - real part of the 1st eigenvalue
692
.  ai     - imaginary part of the 1st eigenvalue
693
.  br     - real part of the 2nd eigenvalue
694
-  bi     - imaginary part of the 2nd eigenvalue
695
 
696
   Output Parameter:
697
.  res    - result of comparison
698
 
699
   Notes:
700
   Returns an integer less than, equal to, or greater than zero if the first
701
   eigenvalue is considered to be respectively less than, equal to, or greater
702
   than the second one.
703
 
704
   The criterion of comparison is related to the 'which' parameter set with
705
   QEPSetWhichEigenpairs().
706
 
707
   Level: developer
708
 
709
.seealso: QEPSortEigenvalues(), QEPSetWhichEigenpairs()
710
@*/
711
PetscErrorCode QEPCompareEigenvalues(QEP qep,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result)
712
{
2317 jroman 713
  PetscReal a,b;
1887 jroman 714
 
715
  PetscFunctionBegin;
716
  switch(qep->which) {
717
    case QEP_LARGEST_MAGNITUDE:
718
    case QEP_SMALLEST_MAGNITUDE:
719
      a = SlepcAbsEigenvalue(ar,ai);
720
      b = SlepcAbsEigenvalue(br,bi);
721
      break;
722
    case QEP_LARGEST_REAL:
723
    case QEP_SMALLEST_REAL:
724
      a = PetscRealPart(ar);
725
      b = PetscRealPart(br);
726
      break;
727
    case QEP_LARGEST_IMAGINARY:
728
    case QEP_SMALLEST_IMAGINARY:
729
#if defined(PETSC_USE_COMPLEX)
730
      a = PetscImaginaryPart(ar);
731
      b = PetscImaginaryPart(br);
732
#else
733
      a = PetscAbsReal(ai);
734
      b = PetscAbsReal(bi);
735
#endif
736
      break;
2214 jroman 737
    default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of which");
1887 jroman 738
  }
739
  switch(qep->which) {
740
    case QEP_LARGEST_MAGNITUDE:
741
    case QEP_LARGEST_REAL:
742
    case QEP_LARGEST_IMAGINARY:
743
      if (a<b) *result = -1;
744
      else if (a>b) *result = 1;
745
      else *result = 0;
746
      break;
747
    default:
748
      if (a>b) *result = -1;
749
      else if (a<b) *result = 1;
750
      else *result = 0;
751
  }
752
  PetscFunctionReturn(0);
753
}
754
 
1896 jroman 755
#undef __FUNCT__  
756
#define __FUNCT__ "QEPGetOperationCounters"
757
/*@
758
   QEPGetOperationCounters - Gets the total number of matrix-vector products, dot
759
   products, and linear solve iterations used by the QEP object during the last
760
   QEPSolve() call.
761
 
762
   Not Collective
763
 
764
   Input Parameter:
765
.  qep - quadratic eigensolver context
766
 
767
   Output Parameter:
768
+  matvecs - number of matrix-vector product operations
769
.  dots    - number of dot product operations
770
-  lits    - number of linear iterations
771
 
772
   Notes:
773
   These counters are reset to zero at each successive call to QEPSolve().
774
 
775
   Level: intermediate
776
 
777
@*/
778
PetscErrorCode QEPGetOperationCounters(QEP qep,PetscInt* matvecs,PetscInt* dots,PetscInt* lits)
779
{
780
  PetscErrorCode ierr;
781
 
782
  PetscFunctionBegin;
2213 jroman 783
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1896 jroman 784
  if (matvecs) *matvecs = qep->matvecs;
785
  if (dots) {
786
    ierr = IPGetOperationCounters(qep->ip,dots);CHKERRQ(ierr);
787
  }
788
  if (lits) *lits = qep->linits;
789
  PetscFunctionReturn(0);
790
}