Actual source code: pepsolve.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    PEP routines related to the solution process

 13:    References:

 15:        [1] C. Campos and J.E. Roman, "Parallel iterative refinement in
 16:            polynomial eigenvalue problems", Numer. Linear Algebra Appl.
 17:            23(4):730-745, 2016.
 18: */

 20: #include <slepc/private/pepimpl.h>
 21: #include <slepc/private/bvimpl.h>
 22: #include <petscdraw.h>

 24: static PetscBool  cited = PETSC_FALSE;
 25: static const char citation[] =
 26:   "@Article{slepc-pep-refine,\n"
 27:   "   author = \"C. Campos and J. E. Roman\",\n"
 28:   "   title = \"Parallel iterative refinement in polynomial eigenvalue problems\",\n"
 29:   "   journal = \"Numer. Linear Algebra Appl.\",\n"
 30:   "   volume = \"23\",\n"
 31:   "   number = \"4\",\n"
 32:   "   pages = \"730--745\",\n"
 33:   "   year = \"2016,\"\n"
 34:   "   doi = \"https://doi.org/10.1002/nla.2052\"\n"
 35:   "}\n";

 37: PetscErrorCode PEPComputeVectors(PEP pep)
 38: {
 39:   PetscFunctionBegin;
 40:   PEPCheckSolved(pep,1);
 41:   if (pep->state==PEP_STATE_SOLVED) PetscTryTypeMethod(pep,computevectors);
 42:   pep->state = PEP_STATE_EIGENVECTORS;
 43:   PetscFunctionReturn(PETSC_SUCCESS);
 44: }

 46: PetscErrorCode PEPExtractVectors(PEP pep)
 47: {
 48:   PetscFunctionBegin;
 49:   PEPCheckSolved(pep,1);
 50:   if (pep->state==PEP_STATE_SOLVED) PetscTryTypeMethod(pep,extractvectors);
 51:   PetscFunctionReturn(PETSC_SUCCESS);
 52: }

 54: /*@
 55:    PEPSolve - Solves the polynomial eigensystem.

 57:    Collective

 59:    Input Parameter:
 60: .  pep - the polynomial eigensolver context

 62:    Options Database Keys:
 63: +  -pep_view - print information about the solver used
 64: .  -pep_view_matk - view the coefficient matrix Ak (replace k by an integer from 0 to nmat-1)
 65: .  -pep_view_vectors - view the computed eigenvectors
 66: .  -pep_view_values - view the computed eigenvalues
 67: .  -pep_converged_reason - print reason for convergence, and number of iterations
 68: .  -pep_error_absolute - print absolute errors of each eigenpair
 69: .  -pep_error_relative - print relative errors of each eigenpair
 70: -  -pep_error_backward - print backward errors of each eigenpair

 72:    Notes:
 73:    All the command-line options listed above admit an optional argument specifying
 74:    the viewer type and options. For instance, use '-pep_view_mat0 binary:amatrix.bin'
 75:    to save the A matrix to a binary file, '-pep_view_values draw' to draw the computed
 76:    eigenvalues graphically, or '-pep_error_relative :myerr.m:ascii_matlab' to save
 77:    the errors in a file that can be executed in Matlab.

 79:    Level: beginner

 81: .seealso: [](ch:pep), `PEPCreate()`, `PEPSetUp()`, `PEPDestroy()`, `PEPSetTolerances()`
 82: @*/
 83: PetscErrorCode PEPSolve(PEP pep)
 84: {
 85:   PetscInt       i,k;
 86:   PetscBool      flg,islinear;
 87:   char           str[16];

 89:   PetscFunctionBegin;
 91:   if (pep->state>=PEP_STATE_SOLVED) PetscFunctionReturn(PETSC_SUCCESS);
 92:   PetscCall(PetscLogEventBegin(PEP_Solve,pep,0,0,0));

 94:   /* call setup */
 95:   PetscCall(PEPSetUp(pep));
 96:   pep->nconv = 0;
 97:   pep->its   = 0;
 98:   k = pep->lineariz? pep->ncv: pep->ncv*(pep->nmat-1);
 99:   for (i=0;i<k;i++) {
100:     pep->eigr[i]   = 0.0;
101:     pep->eigi[i]   = 0.0;
102:     pep->errest[i] = 0.0;
103:     pep->perm[i]   = i;
104:   }
105:   PetscCall(PEPViewFromOptions(pep,NULL,"-pep_view_pre"));
106:   PetscCall(RGViewFromOptions(pep->rg,NULL,"-rg_view"));

108:   /* Call solver */
109:   PetscUseTypeMethod(pep,solve);
110:   PetscCheck(pep->reason,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
111:   pep->state = PEP_STATE_SOLVED;

113:   /* Only the first nconv columns contain useful information */
114:   PetscCall(BVSetActiveColumns(pep->V,0,pep->nconv));

116:   PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&islinear));
117:   if (!islinear) {
118:     PetscCall(STPostSolve(pep->st));
119:     /* Map eigenvalues back to the original problem */
120:     PetscCall(STGetTransform(pep->st,&flg));
121:     if (flg) PetscTryTypeMethod(pep,backtransform);
122:   }

124: #if !defined(PETSC_USE_COMPLEX)
125:   /* reorder conjugate eigenvalues (positive imaginary first) */
126:   for (i=0;i<pep->nconv-1;i++) {
127:     if (pep->eigi[i] != 0) {
128:       if (pep->eigi[i] < 0) {
129:         pep->eigi[i] = -pep->eigi[i];
130:         pep->eigi[i+1] = -pep->eigi[i+1];
131:         /* the next correction only works with eigenvectors */
132:         PetscCall(PEPComputeVectors(pep));
133:         PetscCall(BVScaleColumn(pep->V,i+1,-1.0));
134:       }
135:       i++;
136:     }
137:   }
138: #endif

140:   if (pep->refine!=PEP_REFINE_NONE) PetscCall(PetscCitationsRegister(citation,&cited));

142:   if (pep->refine==PEP_REFINE_SIMPLE && pep->rits>0 && pep->nconv>0) {
143:     PetscCall(PEPComputeVectors(pep));
144:     PetscCall(PEPNewtonRefinementSimple(pep,&pep->rits,pep->rtol,pep->nconv));
145:   }

147:   /* sort eigenvalues according to pep->which parameter */
148:   PetscCall(SlepcSortEigenvalues(pep->sc,pep->nconv,pep->eigr,pep->eigi,pep->perm));
149:   PetscCall(PetscLogEventEnd(PEP_Solve,pep,0,0,0));

151:   /* various viewers */
152:   PetscCall(PEPViewFromOptions(pep,NULL,"-pep_view"));
153:   PetscCall(PEPConvergedReasonViewFromOptions(pep));
154:   PetscCall(PEPErrorViewFromOptions(pep));
155:   PetscCall(PEPValuesViewFromOptions(pep));
156:   PetscCall(PEPVectorsViewFromOptions(pep));
157:   for (i=0;i<pep->nmat;i++) {
158:     PetscCall(PetscSNPrintf(str,sizeof(str),"-pep_view_mat%" PetscInt_FMT,i));
159:     PetscCall(MatViewFromOptions(pep->A[i],(PetscObject)pep,str));
160:   }

162:   /* Remove the initial subspace */
163:   pep->nini = 0;
164:   PetscFunctionReturn(PETSC_SUCCESS);
165: }

167: /*@
168:    PEPGetIterationNumber - Gets the current iteration number. If the
169:    call to PEPSolve() is complete, then it returns the number of iterations
170:    carried out by the solution method.

172:    Not Collective

174:    Input Parameter:
175: .  pep - the polynomial eigensolver context

177:    Output Parameter:
178: .  its - number of iterations

180:    Note:
181:    During the i-th iteration this call returns i-1. If PEPSolve() is
182:    complete, then parameter "its" contains either the iteration number at
183:    which convergence was successfully reached, or failure was detected.
184:    Call PEPGetConvergedReason() to determine if the solver converged or
185:    failed and why.

187:    Level: intermediate

189: .seealso: [](ch:pep), `PEPGetConvergedReason()`, `PEPSetTolerances()`
190: @*/
191: PetscErrorCode PEPGetIterationNumber(PEP pep,PetscInt *its)
192: {
193:   PetscFunctionBegin;
195:   PetscAssertPointer(its,2);
196:   *its = pep->its;
197:   PetscFunctionReturn(PETSC_SUCCESS);
198: }

200: /*@
201:    PEPGetConverged - Gets the number of converged eigenpairs.

203:    Not Collective

205:    Input Parameter:
206: .  pep - the polynomial eigensolver context

208:    Output Parameter:
209: .  nconv - number of converged eigenpairs

211:    Note:
212:    This function should be called after PEPSolve() has finished.

214:    Level: beginner

216: .seealso: [](ch:pep), `PEPSetDimensions()`, `PEPSolve()`, `PEPGetEigenpair()`
217: @*/
218: PetscErrorCode PEPGetConverged(PEP pep,PetscInt *nconv)
219: {
220:   PetscFunctionBegin;
222:   PetscAssertPointer(nconv,2);
223:   PEPCheckSolved(pep,1);
224:   *nconv = pep->nconv;
225:   PetscFunctionReturn(PETSC_SUCCESS);
226: }

228: /*@
229:    PEPGetConvergedReason - Gets the reason why the `PEPSolve()` iteration was
230:    stopped.

232:    Not Collective

234:    Input Parameter:
235: .  pep - the polynomial eigensolver context

237:    Output Parameter:
238: .  reason - negative value indicates diverged, positive value converged, see
239:    `PEPConvergedReason` for the possible values

241:    Options Database Key:
242: .  -pep_converged_reason - print reason for convergence/divergence, and number of iterations

244:    Note:
245:    If this routine is called before or doing the `PEPSolve()` the value of
246:    `PEP_CONVERGED_ITERATING` is returned.

248:    Level: intermediate

250: .seealso: [](ch:pep), `PEPSetTolerances()`, `PEPSolve()`, `PEPConvergedReason`
251: @*/
252: PetscErrorCode PEPGetConvergedReason(PEP pep,PEPConvergedReason *reason)
253: {
254:   PetscFunctionBegin;
256:   PetscAssertPointer(reason,2);
257:   PEPCheckSolved(pep,1);
258:   *reason = pep->reason;
259:   PetscFunctionReturn(PETSC_SUCCESS);
260: }

262: /*@
263:    PEPGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
264:    PEPSolve(). The solution consists in both the eigenvalue and the eigenvector.

266:    Collective

268:    Input Parameters:
269: +  pep - the polynomial eigensolver context
270: -  i   - index of the solution

272:    Output Parameters:
273: +  eigr - real part of eigenvalue
274: .  eigi - imaginary part of eigenvalue
275: .  Vr   - real part of eigenvector
276: -  Vi   - imaginary part of eigenvector

278:    Notes:
279:    It is allowed to pass NULL for Vr and Vi, if the eigenvector is not
280:    required. Otherwise, the caller must provide valid Vec objects, i.e.,
281:    they must be created by the calling program with e.g. MatCreateVecs().

283:    If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
284:    configured with complex scalars the eigenvalue is stored
285:    directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
286:    set to zero). In any case, the user can pass NULL in Vr or Vi if one of
287:    them is not required.

289:    The index i should be a value between 0 and nconv-1 (see PEPGetConverged()).
290:    Eigenpairs are indexed according to the ordering criterion established
291:    with PEPSetWhichEigenpairs().

293:    Level: beginner

295: .seealso: [](ch:pep), `PEPSolve()`, `PEPGetConverged()`, `PEPSetWhichEigenpairs()`
296: @*/
297: PetscErrorCode PEPGetEigenpair(PEP pep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
298: {
299:   PetscInt       k;

301:   PetscFunctionBegin;
306:   PEPCheckSolved(pep,1);
307:   PetscCheck(i>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
308:   PetscCheck(i<pep->nconv,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see PEPGetConverged()");

310:   PetscCall(PEPComputeVectors(pep));
311:   k = pep->perm[i];

313:   /* eigenvalue */
314: #if defined(PETSC_USE_COMPLEX)
315:   if (eigr) *eigr = pep->eigr[k];
316:   if (eigi) *eigi = 0;
317: #else
318:   if (eigr) *eigr = pep->eigr[k];
319:   if (eigi) *eigi = pep->eigi[k];
320: #endif

322:   /* eigenvector */
323:   PetscCall(BV_GetEigenvector(pep->V,k,pep->eigi[k],Vr,Vi));
324:   PetscFunctionReturn(PETSC_SUCCESS);
325: }

327: /*@
328:    PEPGetErrorEstimate - Returns the error estimate associated to the i-th
329:    computed eigenpair.

331:    Not Collective

333:    Input Parameters:
334: +  pep - the polynomial eigensolver context
335: -  i   - index of eigenpair

337:    Output Parameter:
338: .  errest - the error estimate

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

345:    Level: advanced

347: .seealso: [](ch:pep), `PEPComputeError()`
348: @*/
349: PetscErrorCode PEPGetErrorEstimate(PEP pep,PetscInt i,PetscReal *errest)
350: {
351:   PetscFunctionBegin;
353:   PetscAssertPointer(errest,3);
354:   PEPCheckSolved(pep,1);
355:   PetscCheck(i>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
356:   PetscCheck(i<pep->nconv,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see PEPGetConverged()");
357:   *errest = pep->errest[pep->perm[i]];
358:   PetscFunctionReturn(PETSC_SUCCESS);
359: }

361: /*
362:    PEPComputeResidualNorm_Private - Computes the norm of the residual vector
363:    associated with an eigenpair.

365:    Input Parameters:
366:      kr,ki - eigenvalue
367:      xr,xi - eigenvector
368:      z     - array of 4 work vectors (z[2],z[3] not referenced in complex scalars)
369: */
370: PetscErrorCode PEPComputeResidualNorm_Private(PEP pep,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm)
371: {
372:   Mat            *A=pep->A;
373:   PetscInt       i,nmat=pep->nmat;
374:   PetscScalar    t[20],*vals=t,*ivals=NULL;
375:   Vec            u,w;
376: #if !defined(PETSC_USE_COMPLEX)
377:   Vec            ui,wi;
378:   PetscReal      ni;
379:   PetscBool      imag;
380:   PetscScalar    it[20];
381: #endif

383:   PetscFunctionBegin;
384:   u = z[0]; w = z[1];
385:   PetscCall(VecSet(u,0.0));
386: #if !defined(PETSC_USE_COMPLEX)
387:   ui = z[2]; wi = z[3];
388:   ivals = it;
389: #endif
390:   if (nmat>20) {
391:     PetscCall(PetscMalloc1(nmat,&vals));
392: #if !defined(PETSC_USE_COMPLEX)
393:     PetscCall(PetscMalloc1(nmat,&ivals));
394: #endif
395:   }
396:   PetscCall(PEPEvaluateBasis(pep,kr,ki,vals,ivals));
397: #if !defined(PETSC_USE_COMPLEX)
398:   if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON))
399:     imag = PETSC_FALSE;
400:   else {
401:     imag = PETSC_TRUE;
402:     PetscCall(VecSet(ui,0.0));
403:   }
404: #endif
405:   for (i=0;i<nmat;i++) {
406:     if (vals[i]!=0.0) {
407:       PetscCall(MatMult(A[i],xr,w));
408:       PetscCall(VecAXPY(u,vals[i],w));
409:     }
410: #if !defined(PETSC_USE_COMPLEX)
411:     if (imag) {
412:       if (ivals[i]!=0 || vals[i]!=0) {
413:         PetscCall(MatMult(A[i],xi,wi));
414:         if (vals[i]==0) PetscCall(MatMult(A[i],xr,w));
415:       }
416:       if (ivals[i]!=0) {
417:         PetscCall(VecAXPY(u,-ivals[i],wi));
418:         PetscCall(VecAXPY(ui,ivals[i],w));
419:       }
420:       if (vals[i]!=0) PetscCall(VecAXPY(ui,vals[i],wi));
421:     }
422: #endif
423:   }
424:   PetscCall(VecNorm(u,NORM_2,norm));
425: #if !defined(PETSC_USE_COMPLEX)
426:   if (imag) {
427:     PetscCall(VecNorm(ui,NORM_2,&ni));
428:     *norm = SlepcAbsEigenvalue(*norm,ni);
429:   }
430: #endif
431:   if (nmat>20) {
432:     PetscCall(PetscFree(vals));
433: #if !defined(PETSC_USE_COMPLEX)
434:     PetscCall(PetscFree(ivals));
435: #endif
436:   }
437:   PetscFunctionReturn(PETSC_SUCCESS);
438: }

440: /*@
441:    PEPComputeError - Computes the error (based on the residual norm) associated
442:    with the i-th computed eigenpair.

444:    Collective

446:    Input Parameters:
447: +  pep  - the polynomial eigensolver context
448: .  i    - the solution index
449: -  type - the type of error to compute

451:    Output Parameter:
452: .  error - the error

454:    Notes:
455:    The error can be computed in various ways, all of them based on the residual
456:    norm ||P(l)x||_2 where l is the eigenvalue and x is the eigenvector.
457:    See the users guide for additional details.

459:    Level: beginner

461: .seealso: [](ch:pep), `PEPErrorType`, `PEPSolve()`, `PEPGetErrorEstimate()`
462: @*/
463: PetscErrorCode PEPComputeError(PEP pep,PetscInt i,PEPErrorType type,PetscReal *error)
464: {
465:   Vec            xr,xi,w[4];
466:   PetscScalar    kr,ki;
467:   PetscReal      t,z=0.0;
468:   PetscInt       j;
469:   PetscBool      flg;

471:   PetscFunctionBegin;
475:   PetscAssertPointer(error,4);
476:   PEPCheckSolved(pep,1);

478:   /* allocate work vectors */
479: #if defined(PETSC_USE_COMPLEX)
480:   PetscCall(PEPSetWorkVecs(pep,3));
481:   xi   = NULL;
482:   w[2] = NULL;
483:   w[3] = NULL;
484: #else
485:   PetscCall(PEPSetWorkVecs(pep,6));
486:   xi   = pep->work[3];
487:   w[2] = pep->work[4];
488:   w[3] = pep->work[5];
489: #endif
490:   xr   = pep->work[0];
491:   w[0] = pep->work[1];
492:   w[1] = pep->work[2];

494:   /* compute residual norms */
495:   PetscCall(PEPGetEigenpair(pep,i,&kr,&ki,xr,xi));
496:   PetscCall(PEPComputeResidualNorm_Private(pep,kr,ki,xr,xi,w,error));

498:   /* compute error */
499:   switch (type) {
500:     case PEP_ERROR_ABSOLUTE:
501:       break;
502:     case PEP_ERROR_RELATIVE:
503:       *error /= SlepcAbsEigenvalue(kr,ki);
504:       break;
505:     case PEP_ERROR_BACKWARD:
506:       /* initialization of matrix norms */
507:       if (!pep->nrma[pep->nmat-1]) {
508:         for (j=0;j<pep->nmat;j++) {
509:           PetscCall(MatHasOperation(pep->A[j],MATOP_NORM,&flg));
510:           PetscCheck(flg,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
511:           PetscCall(MatNorm(pep->A[j],NORM_INFINITY,&pep->nrma[j]));
512:         }
513:       }
514:       t = SlepcAbsEigenvalue(kr,ki);
515:       for (j=pep->nmat-1;j>=0;j--) {
516:         z = z*t+pep->nrma[j];
517:       }
518:       *error /= z;
519:       break;
520:     default:
521:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
522:   }
523:   PetscFunctionReturn(PETSC_SUCCESS);
524: }