Actual source code: nepsolve.c

slepc-3.16.3 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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

 13:    References:

 15:        [1] C. Campos and J.E. Roman, "NEP: a module for the parallel solution
 16:            of nonlinear eigenvalue problems in SLEPc", ACM Trans. Math. Soft.
 17:            (to appear), arXiv:1910.11712, 2021.
 18: */

 20: #include <slepc/private/nepimpl.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-nep,\n"
 27:   "   author = \"C. Campos and J. E. Roman\",\n"
 28:   "   title = \"{NEP}: a module for the parallel solution of nonlinear eigenvalue problems in {SLEPc}\",\n"
 29:   "   journal = \"{ACM} Trans. Math. Software\",\n"
 30:   "   volume = \"47\",\n"
 31:   "   number = \"3\",\n"
 32:   "   pages = \"23:1--23:29\",\n"
 33:   "   year = \"2021\",\n"
 34:   "   doi = \"10.1145/3447544\"\n"
 35:   "}\n";

 37: PetscErrorCode NEPComputeVectors(NEP nep)
 38: {

 42:   NEPCheckSolved(nep,1);
 43:   if (nep->state==NEP_STATE_SOLVED && nep->ops->computevectors) {
 44:     (*nep->ops->computevectors)(nep);
 45:   }
 46:   nep->state = NEP_STATE_EIGENVECTORS;
 47:   return(0);
 48: }

 50: /*@
 51:    NEPSolve - Solves the nonlinear eigensystem.

 53:    Collective on nep

 55:    Input Parameter:
 56: .  nep - eigensolver context obtained from NEPCreate()

 58:    Options Database Keys:
 59: +  -nep_view - print information about the solver used
 60: .  -nep_view_vectors binary - save the computed eigenvectors to the default binary viewer
 61: .  -nep_view_values - print computed eigenvalues
 62: .  -nep_converged_reason - print reason for convergence, and number of iterations
 63: .  -nep_error_absolute - print absolute errors of each eigenpair
 64: -  -nep_error_relative - print relative errors of each eigenpair

 66:    Level: beginner

 68: .seealso: NEPCreate(), NEPSetUp(), NEPDestroy(), NEPSetTolerances()
 69: @*/
 70: PetscErrorCode NEPSolve(NEP nep)
 71: {
 73:   PetscInt       i;

 77:   if (nep->state>=NEP_STATE_SOLVED) return(0);
 78:   PetscCitationsRegister(citation,&cited);
 79:   PetscLogEventBegin(NEP_Solve,nep,0,0,0);

 81:   /* call setup */
 82:   NEPSetUp(nep);
 83:   nep->nconv = 0;
 84:   nep->its = 0;
 85:   for (i=0;i<nep->ncv;i++) {
 86:     nep->eigr[i]   = 0.0;
 87:     nep->eigi[i]   = 0.0;
 88:     nep->errest[i] = 0.0;
 89:     nep->perm[i]   = i;
 90:   }
 91:   NEPViewFromOptions(nep,NULL,"-nep_view_pre");
 92:   RGViewFromOptions(nep->rg,NULL,"-rg_view");

 94:   /* call solver */
 95:   (*nep->ops->solve)(nep);
 96:   if (!nep->reason) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
 97:   nep->state = NEP_STATE_SOLVED;

 99:   /* Only the first nconv columns contain useful information */
100:   BVSetActiveColumns(nep->V,0,nep->nconv);
101:   if (nep->twosided) { BVSetActiveColumns(nep->W,0,nep->nconv); }

103:   if (nep->refine==NEP_REFINE_SIMPLE && nep->rits>0 && nep->nconv>0) {
104:     NEPComputeVectors(nep);
105:     NEPNewtonRefinementSimple(nep,&nep->rits,nep->rtol,nep->nconv);
106:     nep->state = NEP_STATE_EIGENVECTORS;
107:   }

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

113:   /* various viewers */
114:   NEPViewFromOptions(nep,NULL,"-nep_view");
115:   NEPConvergedReasonViewFromOptions(nep);
116:   NEPErrorViewFromOptions(nep);
117:   NEPValuesViewFromOptions(nep);
118:   NEPVectorsViewFromOptions(nep);

120:   /* Remove the initial subspace */
121:   nep->nini = 0;

123:   /* Reset resolvent information */
124:   MatDestroy(&nep->resolvent);
125:   return(0);
126: }

128: /*@
129:    NEPProjectOperator - Computes the projection of the nonlinear operator.

131:    Collective on nep

133:    Input Parameters:
134: +  nep - the nonlinear eigensolver context
135: .  j0  - initial index
136: -  j1  - final index

138:    Notes:
139:    This is available for split operator only.

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

148:    Level: developer

150: .seealso: NEPSetSplitOperator()
151: @*/
152: PetscErrorCode NEPProjectOperator(NEP nep,PetscInt j0,PetscInt j1)
153: {
155:   PetscInt       k;
156:   Mat            G;

162:   NEPCheckProblem(nep,1);
163:   NEPCheckSplit(nep,1);
164:   BVSetActiveColumns(nep->V,j0,j1);
165:   for (k=0;k<nep->nt;k++) {
166:     DSGetMat(nep->ds,DSMatExtra[k],&G);
167:     BVMatProject(nep->V,nep->A[k],nep->V,G);
168:     DSRestoreMat(nep->ds,DSMatExtra[k],&G);
169:   }
170:   return(0);
171: }

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

176:    Collective on nep

178:    Input Parameters:
179: +  nep    - the nonlinear eigensolver context
180: .  lambda - scalar argument
181: .  x      - vector to be multiplied against
182: -  v      - workspace vector (used only in the case of split form)

184:    Output Parameters:
185: +  y   - result vector
186: .  A   - (optional) Function matrix, for callback interface only
187: -  B   - (unused) preconditioning matrix

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

196:    Level: developer

198: .seealso: NEPSetSplitOperator(), NEPComputeFunction(), NEPApplyAdjoint()
199: @*/
200: PetscErrorCode NEPApplyFunction(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A,Mat B)
201: {
203:   PetscInt       i;
204:   PetscScalar    alpha;


215:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
216:     VecSet(y,0.0);
217:     for (i=0;i<nep->nt;i++) {
218:       FNEvaluateFunction(nep->f[i],lambda,&alpha);
219:       MatMult(nep->A[i],x,v);
220:       VecAXPY(y,alpha,v);
221:     }
222:   } else {
223:     if (!A) A = nep->function;
224:     NEPComputeFunction(nep,lambda,A,A);
225:     MatMult(A,x,y);
226:   }
227:   return(0);
228: }

230: /*@
231:    NEPApplyAdjoint - Applies the adjoint nonlinear function T(lambda)^* to a given vector.

233:    Collective on nep

235:    Input Parameters:
236: +  nep    - the nonlinear eigensolver context
237: .  lambda - scalar argument
238: .  x      - vector to be multiplied against
239: -  v      - workspace vector (used only in the case of split form)

241:    Output Parameters:
242: +  y   - result vector
243: .  A   - (optional) Function matrix, for callback interface only
244: -  B   - (unused) preconditioning matrix

246:    Level: developer

248: .seealso: NEPSetSplitOperator(), NEPComputeFunction(), NEPApplyFunction()
249: @*/
250: PetscErrorCode NEPApplyAdjoint(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A,Mat B)
251: {
253:   PetscInt       i;
254:   PetscScalar    alpha;
255:   Vec            w;


266:   VecDuplicate(x,&w);
267:   VecCopy(x,w);
268:   VecConjugate(w);
269:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
270:     VecSet(y,0.0);
271:     for (i=0;i<nep->nt;i++) {
272:       FNEvaluateFunction(nep->f[i],lambda,&alpha);
273:       MatMultTranspose(nep->A[i],w,v);
274:       VecAXPY(y,alpha,v);
275:     }
276:   } else {
277:     if (!A) A = nep->function;
278:     NEPComputeFunction(nep,lambda,A,A);
279:     MatMultTranspose(A,w,y);
280:   }
281:   VecDestroy(&w);
282:   VecConjugate(y);
283:   return(0);
284: }

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

289:    Collective on nep

291:    Input Parameters:
292: +  nep    - the nonlinear eigensolver context
293: .  lambda - scalar argument
294: .  x      - vector to be multiplied against
295: -  v      - workspace vector (used only in the case of split form)

297:    Output Parameters:
298: +  y   - result vector
299: -  A   - (optional) Jacobian matrix, for callback interface only

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

308:    Level: developer

310: .seealso: NEPSetSplitOperator(), NEPComputeJacobian()
311: @*/
312: PetscErrorCode NEPApplyJacobian(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A)
313: {
315:   PetscInt       i;
316:   PetscScalar    alpha;


326:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
327:     VecSet(y,0.0);
328:     for (i=0;i<nep->nt;i++) {
329:       FNEvaluateDerivative(nep->f[i],lambda,&alpha);
330:       MatMult(nep->A[i],x,v);
331:       VecAXPY(y,alpha,v);
332:     }
333:   } else {
334:     if (!A) A = nep->jacobian;
335:     NEPComputeJacobian(nep,lambda,A);
336:     MatMult(A,x,y);
337:   }
338:   return(0);
339: }

341: /*@
342:    NEPGetIterationNumber - Gets the current iteration number. If the
343:    call to NEPSolve() is complete, then it returns the number of iterations
344:    carried out by the solution method.

346:    Not Collective

348:    Input Parameter:
349: .  nep - the nonlinear eigensolver context

351:    Output Parameter:
352: .  its - number of iterations

354:    Note:
355:    During the i-th iteration this call returns i-1. If NEPSolve() is
356:    complete, then parameter "its" contains either the iteration number at
357:    which convergence was successfully reached, or failure was detected.
358:    Call NEPGetConvergedReason() to determine if the solver converged or
359:    failed and why.

361:    Level: intermediate

363: .seealso: NEPGetConvergedReason(), NEPSetTolerances()
364: @*/
365: PetscErrorCode NEPGetIterationNumber(NEP nep,PetscInt *its)
366: {
370:   *its = nep->its;
371:   return(0);
372: }

374: /*@
375:    NEPGetConverged - Gets the number of converged eigenpairs.

377:    Not Collective

379:    Input Parameter:
380: .  nep - the nonlinear eigensolver context

382:    Output Parameter:
383: .  nconv - number of converged eigenpairs

385:    Note:
386:    This function should be called after NEPSolve() has finished.

388:    Level: beginner

390: .seealso: NEPSetDimensions(), NEPSolve()
391: @*/
392: PetscErrorCode NEPGetConverged(NEP nep,PetscInt *nconv)
393: {
397:   NEPCheckSolved(nep,1);
398:   *nconv = nep->nconv;
399:   return(0);
400: }

402: /*@
403:    NEPGetConvergedReason - Gets the reason why the NEPSolve() iteration was
404:    stopped.

406:    Not Collective

408:    Input Parameter:
409: .  nep - the nonlinear eigensolver context

411:    Output Parameter:
412: .  reason - negative value indicates diverged, positive value converged

414:    Notes:

416:    Possible values for reason are
417: +  NEP_CONVERGED_TOL - converged up to tolerance
418: .  NEP_CONVERGED_USER - converged due to a user-defined condition
419: .  NEP_DIVERGED_ITS - required more than max_it iterations to reach convergence
420: .  NEP_DIVERGED_BREAKDOWN - generic breakdown in method
421: .  NEP_DIVERGED_LINEAR_SOLVE - inner linear solve failed
422: -  NEP_DIVERGED_SUBSPACE_EXHAUSTED - run out of space for the basis in an
423:    unrestarted solver

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

427:    Level: intermediate

429: .seealso: NEPSetTolerances(), NEPSolve(), NEPConvergedReason
430: @*/
431: PetscErrorCode NEPGetConvergedReason(NEP nep,NEPConvergedReason *reason)
432: {
436:   NEPCheckSolved(nep,1);
437:   *reason = nep->reason;
438:   return(0);
439: }

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

445:    Logically Collective on nep

447:    Input Parameters:
448: +  nep - nonlinear eigensolver context
449: -  i   - index of the solution

451:    Output Parameters:
452: +  eigr - real part of eigenvalue
453: .  eigi - imaginary part of eigenvalue
454: .  Vr   - real part of eigenvector
455: -  Vi   - imaginary part of eigenvector

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

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

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

472:    Level: beginner

474: .seealso: NEPSolve(), NEPGetConverged(), NEPSetWhichEigenpairs(), NEPGetLeftEigenvector()
475: @*/
476: PetscErrorCode NEPGetEigenpair(NEP nep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
477: {
478:   PetscInt       k;

486:   NEPCheckSolved(nep,1);
487:   if (i<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
488:   if (i>=nep->nconv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see NEPGetConverged()");

490:   NEPComputeVectors(nep);
491:   k = nep->perm[i];

493:   /* eigenvalue */
494: #if defined(PETSC_USE_COMPLEX)
495:   if (eigr) *eigr = nep->eigr[k];
496:   if (eigi) *eigi = 0;
497: #else
498:   if (eigr) *eigr = nep->eigr[k];
499:   if (eigi) *eigi = nep->eigi[k];
500: #endif

502:   /* eigenvector */
503:   BV_GetEigenvector(nep->V,k,nep->eigi[k],Vr,Vi);
504:   return(0);
505: }

507: /*@
508:    NEPGetLeftEigenvector - Gets the i-th left eigenvector as computed by NEPSolve().

510:    Logically Collective on nep

512:    Input Parameters:
513: +  nep - eigensolver context
514: -  i   - index of the solution

516:    Output Parameters:
517: +  Wr   - real part of left eigenvector
518: -  Wi   - imaginary part of left eigenvector

520:    Notes:
521:    The caller must provide valid Vec objects, i.e., they must be created
522:    by the calling program with e.g. MatCreateVecs().

524:    If the corresponding eigenvalue is real, then Wi is set to zero. If PETSc is
525:    configured with complex scalars the eigenvector is stored directly in Wr
526:    (Wi is set to zero). In any case, the user can pass NULL in Wr or Wi if one of
527:    them is not required.

529:    The index i should be a value between 0 and nconv-1 (see NEPGetConverged()).
530:    Eigensolutions are indexed according to the ordering criterion established
531:    with NEPSetWhichEigenpairs().

533:    Left eigenvectors are available only if the twosided flag was set, see
534:    NEPSetTwoSided().

536:    Level: intermediate

538: .seealso: NEPGetEigenpair(), NEPGetConverged(), NEPSetWhichEigenpairs(), NEPSetTwoSided()
539: @*/
540: PetscErrorCode NEPGetLeftEigenvector(NEP nep,PetscInt i,Vec Wr,Vec Wi)
541: {
543:   PetscInt       k;

550:   NEPCheckSolved(nep,1);
551:   if (!nep->twosided) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONGSTATE,"Must request left vectors with NEPSetTwoSided");
552:   if (i<0 || i>=nep->nconv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
553:   NEPComputeVectors(nep);
554:   k = nep->perm[i];
555:   BV_GetEigenvector(nep->W,k,nep->eigi[k],Wr,Wi);
556:   return(0);
557: }

559: /*@
560:    NEPGetErrorEstimate - Returns the error estimate associated to the i-th
561:    computed eigenpair.

563:    Not Collective

565:    Input Parameters:
566: +  nep - nonlinear eigensolver context
567: -  i   - index of eigenpair

569:    Output Parameter:
570: .  errest - the error estimate

572:    Notes:
573:    This is the error estimate used internally by the eigensolver. The actual
574:    error bound can be computed with NEPComputeError().

576:    Level: advanced

578: .seealso: NEPComputeError()
579: @*/
580: PetscErrorCode NEPGetErrorEstimate(NEP nep,PetscInt i,PetscReal *errest)
581: {
585:   NEPCheckSolved(nep,1);
586:   if (i<0 || i>=nep->nconv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
587:   *errest = nep->errest[nep->perm[i]];
588:   return(0);
589: }

591: /*
592:    NEPComputeResidualNorm_Private - Computes the norm of the residual vector
593:    associated with an eigenpair.

595:    Input Parameters:
596:      adj    - whether the adjoint T^* must be used instead of T
597:      lambda - eigenvalue
598:      x      - eigenvector
599:      w      - array of work vectors (two vectors in split form, one vector otherwise)
600: */
601: PetscErrorCode NEPComputeResidualNorm_Private(NEP nep,PetscBool adj,PetscScalar lambda,Vec x,Vec *w,PetscReal *norm)
602: {
604:   Vec            y,z=NULL;

607:   y = w[0];
608:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) z = w[1];
609:   if (adj) {
610:     NEPApplyAdjoint(nep,lambda,x,z,y,NULL,NULL);
611:   } else {
612:     NEPApplyFunction(nep,lambda,x,z,y,NULL,NULL);
613:   }
614:   VecNorm(y,NORM_2,norm);
615:   return(0);
616: }

618: /*@
619:    NEPComputeError - Computes the error (based on the residual norm) associated
620:    with the i-th computed eigenpair.

622:    Collective on nep

624:    Input Parameters:
625: +  nep  - the nonlinear eigensolver context
626: .  i    - the solution index
627: -  type - the type of error to compute

629:    Output Parameter:
630: .  error - the error

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

637:    Level: beginner

639: .seealso: NEPErrorType, NEPSolve(), NEPGetErrorEstimate()
640: @*/
641: PetscErrorCode NEPComputeError(NEP nep,PetscInt i,NEPErrorType type,PetscReal *error)
642: {
644:   Vec            xr,xi=NULL;
645:   PetscInt       j,nwork,issplit=0;
646:   PetscScalar    kr,ki,s;
647:   PetscReal      er,z=0.0,errorl,nrm;
648:   PetscBool      flg;

655:   NEPCheckSolved(nep,1);

657:   /* allocate work vectors */
658: #if defined(PETSC_USE_COMPLEX)
659:   nwork = 2;
660: #else
661:   nwork = 3;
662: #endif
663:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
664:     issplit = 1;
665:     nwork++;  /* need an extra work vector for NEPComputeResidualNorm_Private */
666:   }
667:   NEPSetWorkVecs(nep,nwork);
668:   xr = nep->work[issplit+1];
669: #if !defined(PETSC_USE_COMPLEX)
670:   xi = nep->work[issplit+2];
671: #endif

673:   /* compute residual norms */
674:   NEPGetEigenpair(nep,i,&kr,&ki,xr,xi);
675: #if !defined(PETSC_USE_COMPLEX)
676:   if (ki) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Not implemented for complex eigenvalues with real scalars");
677: #endif
678:   NEPComputeResidualNorm_Private(nep,PETSC_FALSE,kr,xr,nep->work,error);
679:   VecNorm(xr,NORM_2,&er);

681:   /* if two-sided, compute left residual norm and take the maximum */
682:   if (nep->twosided) {
683:     NEPGetLeftEigenvector(nep,i,xr,xi);
684:     NEPComputeResidualNorm_Private(nep,PETSC_TRUE,kr,xr,nep->work,&errorl);
685:     *error = PetscMax(*error,errorl);
686:   }

688:   /* compute error */
689:   switch (type) {
690:     case NEP_ERROR_ABSOLUTE:
691:       break;
692:     case NEP_ERROR_RELATIVE:
693:       *error /= PetscAbsScalar(kr)*er;
694:       break;
695:     case NEP_ERROR_BACKWARD:
696:       if (nep->fui!=NEP_USER_INTERFACE_SPLIT) {
697:         NEPComputeFunction(nep,kr,nep->function,nep->function);
698:         MatHasOperation(nep->function,MATOP_NORM,&flg);
699:         if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
700:         MatNorm(nep->function,NORM_INFINITY,&nrm);
701:         *error /= nrm*er;
702:         break;
703:       }
704:       /* initialization of matrix norms */
705:       if (!nep->nrma[0]) {
706:         for (j=0;j<nep->nt;j++) {
707:           MatHasOperation(nep->A[j],MATOP_NORM,&flg);
708:           if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
709:           MatNorm(nep->A[j],NORM_INFINITY,&nep->nrma[j]);
710:         }
711:       }
712:       for (j=0;j<nep->nt;j++) {
713:         FNEvaluateFunction(nep->f[j],kr,&s);
714:         z = z + nep->nrma[j]*PetscAbsScalar(s);
715:       }
716:       *error /= z*er;
717:       break;
718:     default:
719:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
720:   }
721:   return(0);
722: }

724: /*@
725:    NEPComputeFunction - Computes the function matrix T(lambda) that has been
726:    set with NEPSetFunction().

728:    Collective on nep

730:    Input Parameters:
731: +  nep    - the NEP context
732: -  lambda - the scalar argument

734:    Output Parameters:
735: +  A   - Function matrix
736: -  B   - optional preconditioning matrix

738:    Notes:
739:    NEPComputeFunction() is typically used within nonlinear eigensolvers
740:    implementations, so most users would not generally call this routine
741:    themselves.

743:    Level: developer

745: .seealso: NEPSetFunction(), NEPGetFunction()
746: @*/
747: PetscErrorCode NEPComputeFunction(NEP nep,PetscScalar lambda,Mat A,Mat B)
748: {
750:   PetscInt       i;
751:   PetscScalar    alpha;

755:   NEPCheckProblem(nep,1);
756:   switch (nep->fui) {
757:   case NEP_USER_INTERFACE_CALLBACK:
758:     if (!nep->computefunction) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_USER,"Must call NEPSetFunction() first");
759:     PetscLogEventBegin(NEP_FunctionEval,nep,A,B,0);
760:     PetscStackPush("NEP user Function function");
761:     (*nep->computefunction)(nep,lambda,A,B,nep->functionctx);
762:     PetscStackPop;
763:     PetscLogEventEnd(NEP_FunctionEval,nep,A,B,0);
764:     break;
765:   case NEP_USER_INTERFACE_SPLIT:
766:     MatZeroEntries(A);
767:     for (i=0;i<nep->nt;i++) {
768:       FNEvaluateFunction(nep->f[i],lambda,&alpha);
769:       MatAXPY(A,alpha,nep->A[i],nep->mstr);
770:     }
771:     if (A != B) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Not implemented");
772:     break;
773:   }
774:   return(0);
775: }

777: /*@
778:    NEPComputeJacobian - Computes the Jacobian matrix T'(lambda) that has been
779:    set with NEPSetJacobian().

781:    Collective on nep

783:    Input Parameters:
784: +  nep    - the NEP context
785: -  lambda - the scalar argument

787:    Output Parameters:
788: .  A   - Jacobian matrix

790:    Notes:
791:    Most users should not need to explicitly call this routine, as it
792:    is used internally within the nonlinear eigensolvers.

794:    Level: developer

796: .seealso: NEPSetJacobian(), NEPGetJacobian()
797: @*/
798: PetscErrorCode NEPComputeJacobian(NEP nep,PetscScalar lambda,Mat A)
799: {
801:   PetscInt       i;
802:   PetscScalar    alpha;

806:   NEPCheckProblem(nep,1);
807:   switch (nep->fui) {
808:   case NEP_USER_INTERFACE_CALLBACK:
809:     if (!nep->computejacobian) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_USER,"Must call NEPSetJacobian() first");
810:     PetscLogEventBegin(NEP_JacobianEval,nep,A,0,0);
811:     PetscStackPush("NEP user Jacobian function");
812:     (*nep->computejacobian)(nep,lambda,A,nep->jacobianctx);
813:     PetscStackPop;
814:     PetscLogEventEnd(NEP_JacobianEval,nep,A,0,0);
815:     break;
816:   case NEP_USER_INTERFACE_SPLIT:
817:     MatZeroEntries(A);
818:     for (i=0;i<nep->nt;i++) {
819:       FNEvaluateDerivative(nep->f[i],lambda,&alpha);
820:       MatAXPY(A,alpha,nep->A[i],nep->mstr);
821:     }
822:     break;
823:   }
824:   return(0);
825: }