Actual source code: nepsolve.c

slepc-3.9.1 2018-05-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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:    NEP routines related to the solution process
 12: */

 14: #include <slepc/private/nepimpl.h>       /*I "slepcnep.h" I*/
 15: #include <petscdraw.h>

 17: PetscErrorCode NEPComputeVectors(NEP nep)
 18: {

 22:   NEPCheckSolved(nep,1);
 23:   if (nep->state==NEP_STATE_SOLVED && nep->ops->computevectors) {
 24:     (*nep->ops->computevectors)(nep);
 25:   }
 26:   nep->state = NEP_STATE_EIGENVECTORS;
 27:   return(0);
 28: }

 30: /*@
 31:    NEPSolve - Solves the nonlinear eigensystem.

 33:    Collective on NEP

 35:    Input Parameter:
 36: .  nep - eigensolver context obtained from NEPCreate()

 38:    Options Database Keys:
 39: +  -nep_view - print information about the solver used
 40: .  -nep_view_vectors binary - save the computed eigenvectors to the default binary viewer
 41: .  -nep_view_values - print computed eigenvalues
 42: .  -nep_converged_reason - print reason for convergence, and number of iterations
 43: .  -nep_error_absolute - print absolute errors of each eigenpair
 44: -  -nep_error_relative - print relative errors of each eigenpair

 46:    Level: beginner

 48: .seealso: NEPCreate(), NEPSetUp(), NEPDestroy(), NEPSetTolerances()
 49: @*/
 50: PetscErrorCode NEPSolve(NEP nep)
 51: {
 53:   PetscInt       i;

 57:   if (nep->state>=NEP_STATE_SOLVED) return(0);
 58:   PetscLogEventBegin(NEP_Solve,nep,0,0,0);

 60:   /* call setup */
 61:   NEPSetUp(nep);
 62:   nep->nconv = 0;
 63:   nep->its = 0;
 64:   for (i=0;i<nep->ncv;i++) {
 65:     nep->eigr[i]   = 0.0;
 66:     nep->eigi[i]   = 0.0;
 67:     nep->errest[i] = 0.0;
 68:     nep->perm[i]   = i;
 69:   }
 70:   NEPViewFromOptions(nep,NULL,"-nep_view_pre");
 71:   RGViewFromOptions(nep->rg,NULL,"-rg_view");

 73:   (*nep->ops->solve)(nep);
 74:   nep->state = NEP_STATE_SOLVED;

 76:   if (!nep->reason) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");

 78:   if (nep->refine==NEP_REFINE_SIMPLE && nep->rits>0 && nep->nconv>0) {
 79:     NEPComputeVectors(nep);
 80:     NEPNewtonRefinementSimple(nep,&nep->rits,nep->rtol,nep->nconv);
 81:     nep->state = NEP_STATE_EIGENVECTORS;
 82:   }

 84:   /* sort eigenvalues according to nep->which parameter */
 85:   SlepcSortEigenvalues(nep->sc,nep->nconv,nep->eigr,nep->eigi,nep->perm);
 86:   PetscLogEventEnd(NEP_Solve,nep,0,0,0);

 88:   /* various viewers */
 89:   NEPViewFromOptions(nep,NULL,"-nep_view");
 90:   NEPReasonViewFromOptions(nep);
 91:   NEPErrorViewFromOptions(nep);
 92:   NEPValuesViewFromOptions(nep);
 93:   NEPVectorsViewFromOptions(nep);

 95:   /* Remove the initial subspace */
 96:   nep->nini = 0;
 97:   return(0);
 98: }

100: /*@
101:    NEPProjectOperator - Computes the projection of the nonlinear operator.

103:    Collective on NEP

105:    Input Parameters:
106: +  nep - the nonlinear eigensolver context
107: .  j0  - initial index
108: -  j1  - final index

110:    Notes:
111:    This is available for split operator only.

113:    The nonlinear operator T(lambda) is projected onto span(V), where V is
114:    an orthonormal basis built internally by the solver. The projected
115:    operator is equal to sum_i V'*A_i*V*f_i(lambda), so this function
116:    computes all matrices Ei = V'*A_i*V, and stores them in the extra
117:    matrices inside DS. Only rows/columns in the range [j0,j1-1] are computed,
118:    the previous ones are assumed to be available already.

120:    Level: developer

122: .seealso: NEPSetSplitOperator()
123: @*/
124: PetscErrorCode NEPProjectOperator(NEP nep,PetscInt j0,PetscInt j1)
125: {
127:   PetscInt       k;
128:   Mat            G;

134:   NEPCheckProblem(nep,1);
135:   NEPCheckSplit(nep,1);
136:   BVSetActiveColumns(nep->V,j0,j1);
137:   for (k=0;k<nep->nt;k++) {
138:     DSGetMat(nep->ds,DSMatExtra[k],&G);
139:     BVMatProject(nep->V,nep->A[k],nep->V,G);
140:     DSRestoreMat(nep->ds,DSMatExtra[k],&G);
141:   }
142:   return(0);
143: }

145: /*@
146:    NEPApplyFunction - Applies the nonlinear function T(lambda) to a given vector.

148:    Collective on NEP

150:    Input Parameters:
151: +  nep    - the nonlinear eigensolver context
152: .  lambda - scalar argument
153: .  x      - vector to be multiplied against
154: -  v      - workspace vector (used only in the case of split form)

156:    Output Parameters:
157: +  y   - result vector
158: .  A   - Function matrix
159: -  B   - optional preconditioning matrix

161:    Note:
162:    If the nonlinear operator is represented in split form, the result
163:    y = T(lambda)*x is computed without building T(lambda) explicitly. In
164:    that case, parameters A and B are not used. Otherwise, the matrix
165:    T(lambda) is built and the effect is the same as a call to
166:    NEPComputeFunction() followed by a MatMult().

168:    Level: developer

170: .seealso: NEPSetSplitOperator(), NEPComputeFunction()
171: @*/
172: PetscErrorCode NEPApplyFunction(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A,Mat B)
173: {
175:   PetscInt       i;
176:   PetscScalar    alpha;


187:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
188:     VecSet(y,0.0);
189:     for (i=0;i<nep->nt;i++) {
190:       FNEvaluateFunction(nep->f[i],lambda,&alpha);
191:       MatMult(nep->A[i],x,v);
192:       VecAXPY(y,alpha,v);
193:     }
194:   } else {
195:     NEPComputeFunction(nep,lambda,A,B);
196:     MatMult(A,x,y);
197:   }
198:   return(0);
199: }

201: /*@
202:    NEPApplyJacobian - Applies the nonlinear Jacobian T'(lambda) to a given vector.

204:    Collective on NEP

206:    Input Parameters:
207: +  nep    - the nonlinear eigensolver context
208: .  lambda - scalar argument
209: .  x      - vector to be multiplied against
210: -  v      - workspace vector (used only in the case of split form)

212:    Output Parameters:
213: +  y   - result vector
214: -  A   - Jacobian matrix

216:    Note:
217:    If the nonlinear operator is represented in split form, the result
218:    y = T'(lambda)*x is computed without building T'(lambda) explicitly. In
219:    that case, parameter A is not used. Otherwise, the matrix
220:    T'(lambda) is built and the effect is the same as a call to
221:    NEPComputeJacobian() followed by a MatMult().

223:    Level: developer

225: .seealso: NEPSetSplitOperator(), NEPComputeJacobian()
226: @*/
227: PetscErrorCode NEPApplyJacobian(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A)
228: {
230:   PetscInt       i;
231:   PetscScalar    alpha;


241:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
242:     VecSet(y,0.0);
243:     for (i=0;i<nep->nt;i++) {
244:       FNEvaluateDerivative(nep->f[i],lambda,&alpha);
245:       MatMult(nep->A[i],x,v);
246:       VecAXPY(y,alpha,v);
247:     }
248:   } else {
249:     NEPComputeJacobian(nep,lambda,A);
250:     MatMult(A,x,y);
251:   }
252:   return(0);
253: }

255: /*@
256:    NEPGetIterationNumber - Gets the current iteration number. If the
257:    call to NEPSolve() is complete, then it returns the number of iterations
258:    carried out by the solution method.

260:    Not Collective

262:    Input Parameter:
263: .  nep - the nonlinear eigensolver context

265:    Output Parameter:
266: .  its - number of iterations

268:    Note:
269:    During the i-th iteration this call returns i-1. If NEPSolve() is
270:    complete, then parameter "its" contains either the iteration number at
271:    which convergence was successfully reached, or failure was detected.
272:    Call NEPGetConvergedReason() to determine if the solver converged or
273:    failed and why.

275:    Level: intermediate

277: .seealso: NEPGetConvergedReason(), NEPSetTolerances()
278: @*/
279: PetscErrorCode NEPGetIterationNumber(NEP nep,PetscInt *its)
280: {
284:   *its = nep->its;
285:   return(0);
286: }

288: /*@
289:    NEPGetConverged - Gets the number of converged eigenpairs.

291:    Not Collective

293:    Input Parameter:
294: .  nep - the nonlinear eigensolver context

296:    Output Parameter:
297: .  nconv - number of converged eigenpairs

299:    Note:
300:    This function should be called after NEPSolve() has finished.

302:    Level: beginner

304: .seealso: NEPSetDimensions(), NEPSolve()
305: @*/
306: PetscErrorCode NEPGetConverged(NEP nep,PetscInt *nconv)
307: {
311:   NEPCheckSolved(nep,1);
312:   *nconv = nep->nconv;
313:   return(0);
314: }

316: /*@
317:    NEPGetConvergedReason - Gets the reason why the NEPSolve() iteration was
318:    stopped.

320:    Not Collective

322:    Input Parameter:
323: .  nep - the nonlinear eigensolver context

325:    Output Parameter:
326: .  reason - negative value indicates diverged, positive value converged

328:    Notes:

330:    Possible values for reason are
331: +  NEP_CONVERGED_TOL - converged up to tolerance
332: .  NEP_CONVERGED_USER - converged due to a user-defined condition
333: .  NEP_DIVERGED_ITS - required more than max_it iterations to reach convergence
334: .  NEP_DIVERGED_BREAKDOWN - generic breakdown in method
335: -  NEP_DIVERGED_LINEAR_SOLVE - inner linear solve failed

337:    Can only be called after the call to NEPSolve() is complete.

339:    Level: intermediate

341: .seealso: NEPSetTolerances(), NEPSolve(), NEPConvergedReason
342: @*/
343: PetscErrorCode NEPGetConvergedReason(NEP nep,NEPConvergedReason *reason)
344: {
348:   NEPCheckSolved(nep,1);
349:   *reason = nep->reason;
350:   return(0);
351: }

353: /*@C
354:    NEPGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
355:    NEPSolve(). The solution consists in both the eigenvalue and the eigenvector.

357:    Logically Collective on NEP

359:    Input Parameters:
360: +  nep - nonlinear eigensolver context
361: -  i   - index of the solution

363:    Output Parameters:
364: +  eigr - real part of eigenvalue
365: .  eigi - imaginary part of eigenvalue
366: .  Vr   - real part of eigenvector
367: -  Vi   - imaginary part of eigenvector

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

374:    If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
375:    configured with complex scalars the eigenvalue is stored
376:    directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
377:    set to zero). In both cases, the user can pass NULL in eigi and Vi.

379:    The index i should be a value between 0 and nconv-1 (see NEPGetConverged()).
380:    Eigenpairs are indexed according to the ordering criterion established
381:    with NEPSetWhichEigenpairs().

383:    Level: beginner

385: .seealso: NEPSolve(), NEPGetConverged(), NEPSetWhichEigenpairs()
386: @*/
387: PetscErrorCode NEPGetEigenpair(NEP nep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
388: {
389:   PetscInt       k;

397:   NEPCheckSolved(nep,1);
398:   if (i<0 || i>=nep->nconv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");

400:   NEPComputeVectors(nep);
401:   k = nep->perm[i];

403:   /* eigenvalue */
404: #if defined(PETSC_USE_COMPLEX)
405:   if (eigr) *eigr = nep->eigr[k];
406:   if (eigi) *eigi = 0;
407: #else
408:   if (eigr) *eigr = nep->eigr[k];
409:   if (eigi) *eigi = nep->eigi[k];
410: #endif

412:   /* eigenvector */
413: #if defined(PETSC_USE_COMPLEX)
414:   if (Vr) { BVCopyVec(nep->V,k,Vr); }
415:   if (Vi) { VecSet(Vi,0.0); }
416: #else
417:   if (nep->eigi[k]>0) { /* first value of conjugate pair */
418:     if (Vr) { BVCopyVec(nep->V,k,Vr); }
419:     if (Vi) { BVCopyVec(nep->V,k+1,Vi); }
420:   } else if (nep->eigi[k]<0) { /* second value of conjugate pair */
421:     if (Vr) { BVCopyVec(nep->V,k-1,Vr); }
422:     if (Vi) {
423:       BVCopyVec(nep->V,k,Vi);
424:       VecScale(Vi,-1.0);
425:     }
426:   } else { /* real eigenvalue */
427:     if (Vr) { BVCopyVec(nep->V,k,Vr); }
428:     if (Vi) { VecSet(Vi,0.0); }
429:   }
430: #endif
431:   return(0);
432: }

434: /*@
435:    NEPGetErrorEstimate - Returns the error estimate associated to the i-th
436:    computed eigenpair.

438:    Not Collective

440:    Input Parameter:
441: +  nep - nonlinear eigensolver context
442: -  i   - index of eigenpair

444:    Output Parameter:
445: .  errest - the error estimate

447:    Notes:
448:    This is the error estimate used internally by the eigensolver. The actual
449:    error bound can be computed with NEPComputeRelativeError().

451:    Level: advanced

453: .seealso: NEPComputeRelativeError()
454: @*/
455: PetscErrorCode NEPGetErrorEstimate(NEP nep,PetscInt i,PetscReal *errest)
456: {
460:   NEPCheckSolved(nep,1);
461:   if (i<0 || i>=nep->nconv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
462:   *errest = nep->errest[nep->perm[i]];
463:   return(0);
464: }

466: /*
467:    NEPComputeResidualNorm_Private - Computes the norm of the residual vector
468:    associated with an eigenpair.

470:    Input Parameters:
471:      lambda - eigenvalue
472:      x      - eigenvector
473:      w      - array of work vectors (two vectors in split form, one vector otherwise)
474: */
475: PetscErrorCode NEPComputeResidualNorm_Private(NEP nep,PetscScalar lambda,Vec x,Vec *w,PetscReal *norm)
476: {
478:   Vec            y,z=NULL;

481:   y = w[0];
482:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) z = w[1];
483:   NEPApplyFunction(nep,lambda,x,z,y,nep->function,nep->function_pre);
484:   VecNorm(y,NORM_2,norm);
485:   return(0);
486: }

488: /*@
489:    NEPComputeError - Computes the error (based on the residual norm) associated
490:    with the i-th computed eigenpair.

492:    Collective on NEP

494:    Input Parameter:
495: +  nep  - the nonlinear eigensolver context
496: .  i    - the solution index
497: -  type - the type of error to compute

499:    Output Parameter:
500: .  error - the error

502:    Notes:
503:    The error can be computed in various ways, all of them based on the residual
504:    norm computed as ||T(lambda)x||_2 where lambda is the eigenvalue and x is the
505:    eigenvector.

507:    Level: beginner

509: .seealso: NEPErrorType, NEPSolve(), NEPGetErrorEstimate()
510: @*/
511: PetscErrorCode NEPComputeError(NEP nep,PetscInt i,NEPErrorType type,PetscReal *error)
512: {
514:   Vec            xr,xi=NULL;
515:   PetscInt       j,nwork,issplit=0;
516:   PetscScalar    kr,ki,s;
517:   PetscReal      er,z=0.0;
518:   PetscBool      flg;

525:   NEPCheckSolved(nep,1);

527:   /* allocate work vectors */
528: #if defined(PETSC_USE_COMPLEX)
529:   nwork = 2;
530: #else
531:   nwork = 3;
532: #endif
533:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
534:     issplit = 1;
535:     nwork++;  /* need an extra work vector for NEPComputeResidualNorm_Private */
536:   }
537:   NEPSetWorkVecs(nep,nwork);
538:   xr = nep->work[issplit+1];
539: #if !defined(PETSC_USE_COMPLEX)
540:   xi = nep->work[issplit+2];
541: #endif

543:   /* compute residual norms */
544:   NEPGetEigenpair(nep,i,&kr,&ki,xr,xi);
545: #if !defined(PETSC_USE_COMPLEX)
546:   if (ki) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Not implemented for complex eigenvalues with real scalars");
547: #endif
548:   NEPComputeResidualNorm_Private(nep,kr,xr,nep->work,error);
549:   VecNorm(xr,NORM_2,&er);

551:   /* compute error */
552:   switch (type) {
553:     case NEP_ERROR_ABSOLUTE:
554:       break;
555:     case NEP_ERROR_RELATIVE:
556:       *error /= PetscAbsScalar(kr)*er;
557:       break;
558:     case NEP_ERROR_BACKWARD:
559:       if (nep->fui!=NEP_USER_INTERFACE_SPLIT) {
560:         *error = 0.0;
561:         PetscInfo(nep,"Backward error only available in split form\n");
562:         break;
563:       }
564:       /* initialization of matrix norms */
565:       if (!nep->nrma[0]) {
566:         for (j=0;j<nep->nt;j++) {
567:           MatHasOperation(nep->A[j],MATOP_NORM,&flg);
568:           if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
569:           MatNorm(nep->A[j],NORM_INFINITY,&nep->nrma[j]);
570:         }
571:       }
572:       for (j=0;j<nep->nt;j++) {
573:         FNEvaluateFunction(nep->f[j],kr,&s);
574:         z = z + nep->nrma[j]*PetscAbsScalar(s);
575:       }
576:       *error /= z;
577:       break;
578:     default:
579:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
580:   }
581:   return(0);
582: }

584: /*@
585:    NEPComputeFunction - Computes the function matrix T(lambda) that has been
586:    set with NEPSetFunction().

588:    Collective on NEP and Mat

590:    Input Parameters:
591: +  nep    - the NEP context
592: -  lambda - the scalar argument

594:    Output Parameters:
595: +  A   - Function matrix
596: -  B   - optional preconditioning matrix

598:    Notes:
599:    NEPComputeFunction() is typically used within nonlinear eigensolvers
600:    implementations, so most users would not generally call this routine
601:    themselves.

603:    Level: developer

605: .seealso: NEPSetFunction(), NEPGetFunction()
606: @*/
607: PetscErrorCode NEPComputeFunction(NEP nep,PetscScalar lambda,Mat A,Mat B)
608: {
610:   PetscInt       i;
611:   PetscScalar    alpha;

615:   NEPCheckProblem(nep,1);
616:   switch (nep->fui) {
617:   case NEP_USER_INTERFACE_CALLBACK:
618:     if (!nep->computefunction) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_USER,"Must call NEPSetFunction() first");
619:     PetscLogEventBegin(NEP_FunctionEval,nep,A,B,0);
620:     PetscStackPush("NEP user Function function");
621:     (*nep->computefunction)(nep,lambda,A,B,nep->functionctx);
622:     PetscStackPop;
623:     PetscLogEventEnd(NEP_FunctionEval,nep,A,B,0);
624:     break;
625:   case NEP_USER_INTERFACE_SPLIT:
626:     MatZeroEntries(A);
627:     for (i=0;i<nep->nt;i++) {
628:       FNEvaluateFunction(nep->f[i],lambda,&alpha);
629:       MatAXPY(A,alpha,nep->A[i],nep->mstr);
630:     }
631:     if (A != B) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Not implemented");
632:     break;
633:   case NEP_USER_INTERFACE_DERIVATIVES:
634:     PetscLogEventBegin(NEP_DerivativesEval,nep,A,B,0);
635:     PetscStackPush("NEP user Derivatives function");
636:     (*nep->computederivatives)(nep,lambda,0,A,nep->derivativesctx);
637:     PetscStackPop;
638:     PetscLogEventEnd(NEP_DerivativesEval,nep,A,B,0);
639:     break;
640:   }
641:   return(0);
642: }

644: /*@
645:    NEPComputeJacobian - Computes the Jacobian matrix T'(lambda) that has been
646:    set with NEPSetJacobian().

648:    Collective on NEP and Mat

650:    Input Parameters:
651: +  nep    - the NEP context
652: -  lambda - the scalar argument

654:    Output Parameters:
655: .  A   - Jacobian matrix

657:    Notes:
658:    Most users should not need to explicitly call this routine, as it
659:    is used internally within the nonlinear eigensolvers.

661:    Level: developer

663: .seealso: NEPSetJacobian(), NEPGetJacobian()
664: @*/
665: PetscErrorCode NEPComputeJacobian(NEP nep,PetscScalar lambda,Mat A)
666: {
668:   PetscInt       i;
669:   PetscScalar    alpha;

673:   NEPCheckProblem(nep,1);
674:   switch (nep->fui) {
675:   case NEP_USER_INTERFACE_CALLBACK:
676:     if (!nep->computejacobian) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_USER,"Must call NEPSetJacobian() first");
677:     PetscLogEventBegin(NEP_JacobianEval,nep,A,0,0);
678:     PetscStackPush("NEP user Jacobian function");
679:     (*nep->computejacobian)(nep,lambda,A,nep->jacobianctx);
680:     PetscStackPop;
681:     PetscLogEventEnd(NEP_JacobianEval,nep,A,0,0);
682:     break;
683:   case NEP_USER_INTERFACE_SPLIT:
684:     MatZeroEntries(A);
685:     for (i=0;i<nep->nt;i++) {
686:       FNEvaluateDerivative(nep->f[i],lambda,&alpha);
687:       MatAXPY(A,alpha,nep->A[i],nep->mstr);
688:     }
689:     break;
690:   case NEP_USER_INTERFACE_DERIVATIVES:
691:     PetscLogEventBegin(NEP_DerivativesEval,nep,A,0,0);
692:     PetscStackPush("NEP user Derivatives function");
693:     (*nep->computederivatives)(nep,lambda,1,A,nep->derivativesctx);
694:     PetscStackPop;
695:     PetscLogEventEnd(NEP_DerivativesEval,nep,A,0,0);
696:     break;
697:   }
698:   return(0);
699: }