Actual source code: nepsolve.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:       NEP routines related to the solution process.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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.

 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.

 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: */

 24: #include <slepc-private/nepimpl.h>       /*I "slepcnep.h" I*/
 25: #include <petscdraw.h>

 29: PETSC_STATIC_INLINE PetscErrorCode NEPComputeVectors(NEP nep)
 30: {

 34:   NEPCheckSolved(nep,1);
 35:   switch (nep->state) {
 36:   case NEP_STATE_SOLVED:
 37:     if (nep->ops->computevectors) {
 38:       (*nep->ops->computevectors)(nep);
 39:     }
 40:     break;
 41:   default:
 42:     break;
 43:   }
 44:   nep->state = NEP_STATE_EIGENVECTORS;
 45:   return(0);
 46: }

 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_plot_eigs - plot computed eigenvalues

 62:    Level: beginner

 64: .seealso: NEPCreate(), NEPSetUp(), NEPDestroy(), NEPSetTolerances()
 65: @*/
 66: PetscErrorCode NEPSolve(NEP nep)
 67: {
 68:   PetscErrorCode    ierr;
 69:   PetscInt          i;
 70:   PetscReal         re,im;
 71:   PetscBool         flg;
 72:   PetscViewer       viewer;
 73:   PetscViewerFormat format;
 74:   PetscDraw         draw;
 75:   PetscDrawSP       drawsp;

 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:   }
 90:   nep->ktol = 0.1;
 91:   NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->errest,nep->ncv);

 93:   (*nep->ops->solve)(nep);

 95:   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:   if (nep->refine==NEP_REFINE_SIMPLE && nep->rits>0) {
100:     NEPComputeVectors(nep);
101:     NEPNewtonRefinementSimple(nep,&nep->rits,&nep->reftol,nep->nconv);
102:     nep->state = NEP_STATE_EIGENVECTORS;
103:   }

105:   /* sort eigenvalues according to nep->which parameter */
106:   SlepcSortEigenvalues(nep->sc,nep->nconv,nep->eigr,nep->eigi,nep->perm);

108:   PetscLogEventEnd(NEP_Solve,nep,0,0,0);

110:   /* various viewers */
111:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->prefix,"-nep_view",&viewer,&format,&flg);
112:   if (flg && !PetscPreLoadingOn) {
113:     PetscViewerPushFormat(viewer,format);
114:     NEPView(nep,viewer);
115:     PetscViewerPopFormat(viewer);
116:     PetscViewerDestroy(&viewer);
117:   }

119:   flg = PETSC_FALSE;
120:   PetscOptionsGetBool(((PetscObject)nep)->prefix,"-nep_plot_eigs",&flg,NULL);
121:   if (flg) {
122:     PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
123:     PetscViewerDrawGetDraw(viewer,0,&draw);
124:     PetscDrawSPCreate(draw,1,&drawsp);
125:     for (i=0;i<nep->nconv;i++) {
126: #if defined(PETSC_USE_COMPLEX)
127:       re = PetscRealPart(nep->eigr[i]);
128:       im = PetscImaginaryPart(nep->eigi[i]);
129: #else
130:       re = nep->eigr[i];
131:       im = nep->eigi[i];
132: #endif
133:       PetscDrawSPAddPoint(drawsp,&re,&im);
134:     }
135:     PetscDrawSPDraw(drawsp,PETSC_TRUE);
136:     PetscDrawSPDestroy(&drawsp);
137:     PetscViewerDestroy(&viewer);
138:   }

140:   /* Remove the initial subspace */
141:   nep->nini = 0;
142:   return(0);
143: }

147: /*@
148:    NEPProjectOperator - Computes the projection of the nonlinear operator.

150:    Collective on NEP

152:    Input Parameters:
153: +  nep - the nonlinear eigensolver context
154: .  j0  - initial index
155: -  j1  - final index

157:    Notes:
158:    This is available for split operator only.

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

167:    Level: developer

169: .seealso: NEPSetSplitOperator()
170: @*/
171: PetscErrorCode NEPProjectOperator(NEP nep,PetscInt j0,PetscInt j1)
172: {
174:   PetscInt       k;
175:   Mat            G;

181:   if (!nep->split) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONGSTATE,"This solver requires a split operator");
182:   BVSetActiveColumns(nep->V,j0,j1);
183:   for (k=0;k<nep->nt;k++) {
184:     DSGetMat(nep->ds,DSMatExtra[k],&G);
185:     BVMatProject(nep->V,nep->A[k],nep->V,G);
186:     DSRestoreMat(nep->ds,DSMatExtra[k],&G);
187:   }
188:   return(0);
189: }

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

196:    Collective on NEP

198:    Input Parameters:
199: +  nep    - the nonlinear eigensolver context
200: .  lambda - scalar argument
201: .  x      - vector to be multiplied against
202: -  v      - workspace vector

204:    Output Parameters:
205: +  y   - result vector
206: .  A   - Function matrix
207: -  B   - optional preconditioning matrix

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

216:    Level: developer

218: .seealso: NEPSetSplitOperator(), NEPComputeFunction()
219: @*/
220: PetscErrorCode NEPApplyFunction(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A,Mat B)
221: {
223:   PetscInt       i;
224:   PetscScalar    alpha;

232:   if (nep->split) {
233:     VecZeroEntries(y);
234:     for (i=0;i<nep->nt;i++) {
235:       FNEvaluateFunction(nep->f[i],lambda,&alpha);
236:       MatMult(nep->A[i],x,v);
237:       VecAXPY(y,alpha,v);
238:     }
239:   } else {
240:     NEPComputeFunction(nep,lambda,A,B);
241:     MatMult(A,x,y);
242:   }
243:   return(0);
244: }

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

251:    Collective on NEP

253:    Input Parameters:
254: +  nep    - the nonlinear eigensolver context
255: .  lambda - scalar argument
256: .  x      - vector to be multiplied against
257: -  v      - workspace vector

259:    Output Parameters:
260: +  y   - result vector
261: -  A   - Jacobian matrix

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

270:    Level: developer

272: .seealso: NEPSetSplitOperator(), NEPComputeJacobian()
273: @*/
274: PetscErrorCode NEPApplyJacobian(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A)
275: {
277:   PetscInt       i;
278:   PetscScalar    alpha;

286:   if (nep->split) {
287:     VecZeroEntries(y);
288:     for (i=0;i<nep->nt;i++) {
289:       FNEvaluateDerivative(nep->f[i],lambda,&alpha);
290:       MatMult(nep->A[i],x,v);
291:       VecAXPY(y,alpha,v);
292:     }
293:   } else {
294:     NEPComputeJacobian(nep,lambda,A);
295:     MatMult(A,x,y);
296:   }
297:   return(0);
298: }

302: /*@
303:    NEPGetIterationNumber - Gets the current iteration number. If the
304:    call to NEPSolve() is complete, then it returns the number of iterations
305:    carried out by the solution method.

307:    Not Collective

309:    Input Parameter:
310: .  nep - the nonlinear eigensolver context

312:    Output Parameter:
313: .  its - number of iterations

315:    Level: intermediate

317:    Note:
318:    During the i-th iteration this call returns i-1. If NEPSolve() is
319:    complete, then parameter "its" contains either the iteration number at
320:    which convergence was successfully reached, or failure was detected.
321:    Call NEPGetConvergedReason() to determine if the solver converged or
322:    failed and why.

324: .seealso: NEPGetConvergedReason(), NEPSetTolerances()
325: @*/
326: PetscErrorCode NEPGetIterationNumber(NEP nep,PetscInt *its)
327: {
331:   *its = nep->its;
332:   return(0);
333: }

337: /*@
338:    NEPGetConverged - Gets the number of converged eigenpairs.

340:    Not Collective

342:    Input Parameter:
343: .  nep - the nonlinear eigensolver context

345:    Output Parameter:
346: .  nconv - number of converged eigenpairs

348:    Note:
349:    This function should be called after NEPSolve() has finished.

351:    Level: beginner

353: .seealso: NEPSetDimensions(), NEPSolve()
354: @*/
355: PetscErrorCode NEPGetConverged(NEP nep,PetscInt *nconv)
356: {
360:   NEPCheckSolved(nep,1);
361:   *nconv = nep->nconv;
362:   return(0);
363: }

367: /*@C
368:    NEPGetConvergedReason - Gets the reason why the NEPSolve() iteration was
369:    stopped.

371:    Not Collective

373:    Input Parameter:
374: .  nep - the nonlinear eigensolver context

376:    Output Parameter:
377: .  reason - negative value indicates diverged, positive value converged

379:    Possible values for reason:
380: +  NEP_CONVERGED_FNORM_ABS - function norm satisfied absolute tolerance
381: .  NEP_CONVERGED_FNORM_RELATIVE - function norm satisfied relative tolerance
382: .  NEP_CONVERGED_SNORM_RELATIVE - step norm satisfied relative tolerance
383: .  NEP_DIVERGED_LINEAR_SOLVE - inner linear solve failed
384: .  NEP_DIVERGED_FUNCTION_COUNT - reached maximum allowed function evaluations
385: .  NEP_DIVERGED_MAX_IT - required more than its to reach convergence
386: .  NEP_DIVERGED_BREAKDOWN - generic breakdown in method
387: -  NEP_DIVERGED_FNORM_NAN - Inf or NaN detected in function evaluation

389:    Note:
390:    Can only be called after the call to NEPSolve() is complete.

392:    Level: intermediate

394: .seealso: NEPSetTolerances(), NEPSolve(), NEPConvergedReason
395: @*/
396: PetscErrorCode NEPGetConvergedReason(NEP nep,NEPConvergedReason *reason)
397: {
401:   NEPCheckSolved(nep,1);
402:   *reason = nep->reason;
403:   return(0);
404: }

408: /*@
409:    NEPGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
410:    NEPSolve(). The solution consists in both the eigenvalue and the eigenvector.

412:    Logically Collective on NEP

414:    Input Parameters:
415: +  nep - nonlinear eigensolver context
416: -  i   - index of the solution

418:    Output Parameters:
419: +  eigr - real part of eigenvalue
420: .  eigi - imaginary part of eigenvalue
421: .  Vr   - real part of eigenvector
422: -  Vi   - imaginary part of eigenvector

424:    Notes:
425:    If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
426:    configured with complex scalars the eigenvalue is stored
427:    directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
428:    set to zero).

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

434:    Level: beginner

436: .seealso: NEPSolve(), NEPGetConverged(), NEPSetWhichEigenpairs()
437: @*/
438: PetscErrorCode NEPGetEigenpair(NEP nep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
439: {
440:   PetscInt       k;

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

451:   NEPComputeVectors(nep);
452:   if (!nep->perm) k = i;
453:   else k = nep->perm[i];

455:   /* eigenvalue */
456: #if defined(PETSC_USE_COMPLEX)
457:   if (eigr) *eigr = nep->eigr[k];
458:   if (eigi) *eigi = 0;
459: #else
460:   if (eigr) *eigr = nep->eigr[k];
461:   if (eigi) *eigi = nep->eigi[k];
462: #endif

464:   /* eigenvector */
465: #if defined(PETSC_USE_COMPLEX)
466:   if (Vr) { BVCopyVec(nep->V,k,Vr); }
467:   if (Vi) { VecSet(Vi,0.0); }
468: #else
469:   if (nep->eigi[k]>0) { /* first value of conjugate pair */
470:     if (Vr) { BVCopyVec(nep->V,k,Vr); }
471:     if (Vi) { BVCopyVec(nep->V,k+1,Vi); }
472:   } else if (nep->eigi[k]<0) { /* second value of conjugate pair */
473:     if (Vr) { BVCopyVec(nep->V,k-1,Vr); }
474:     if (Vi) {
475:       BVCopyVec(nep->V,k,Vi);
476:       VecScale(Vi,-1.0);
477:     }
478:   } else { /* real eigenvalue */
479:     if (Vr) { BVCopyVec(nep->V,k,Vr); }
480:     if (Vi) { VecSet(Vi,0.0); }
481:   }
482: #endif
483:   return(0);
484: }

488: /*@
489:    NEPGetErrorEstimate - Returns the error estimate associated to the i-th
490:    computed eigenpair.

492:    Not Collective

494:    Input Parameter:
495: +  nep - nonlinear eigensolver context
496: -  i   - index of eigenpair

498:    Output Parameter:
499: .  errest - the error estimate

501:    Notes:
502:    This is the error estimate used internally by the eigensolver. The actual
503:    error bound can be computed with NEPComputeRelativeError().

505:    Level: advanced

507: .seealso: NEPComputeRelativeError()
508: @*/
509: PetscErrorCode NEPGetErrorEstimate(NEP nep,PetscInt i,PetscReal *errest)
510: {
514:   NEPCheckSolved(nep,1);
515:   if (i<0 || i>=nep->nconv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
516:   if (nep->perm) i = nep->perm[i];
517:   if (errest) *errest = nep->errest[i];
518:   return(0);
519: }

523: /*
524:    NEPComputeResidualNorm_Private - Computes the norm of the residual vector
525:    associated with an eigenpair.
526: */
527: PetscErrorCode NEPComputeResidualNorm_Private(NEP nep,PetscScalar lambda,Vec x,PetscReal *norm)
528: {
530:   Vec            u;
531:   Mat            T=nep->function;

534:   BVGetVec(nep->V,&u);
535:   NEPComputeFunction(nep,lambda,T,T);
536:   MatMult(T,x,u);
537:   VecNorm(u,NORM_2,norm);
538:   VecDestroy(&u);
539:   return(0);
540: }

544: /*@
545:    NEPComputeResidualNorm - Computes the norm of the residual vector associated with
546:    the i-th computed eigenpair.

548:    Collective on NEP

550:    Input Parameter:
551: +  nep - the nonlinear eigensolver context
552: -  i   - the solution index

554:    Output Parameter:
555: .  norm - the residual norm, computed as ||T(lambda)x||_2 where lambda is the
556:    eigenvalue and x is the eigenvector.

558:    Notes:
559:    The index i should be a value between 0 and nconv-1 (see NEPGetConverged()).
560:    Eigenpairs are indexed according to the ordering criterion established
561:    with NEPSetWhichEigenpairs().

563:    Level: beginner

565: .seealso: NEPSolve(), NEPGetConverged(), NEPSetWhichEigenpairs()
566: @*/
567: PetscErrorCode NEPComputeResidualNorm(NEP nep,PetscInt i,PetscReal *norm)
568: {
570:   Vec            xr,xi;
571:   PetscScalar    kr,ki;

577:   NEPCheckSolved(nep,1);
578:   BVGetVec(nep->V,&xr);
579:   BVGetVec(nep->V,&xi);
580:   NEPGetEigenpair(nep,i,&kr,&ki,xr,xi);
581: #if !defined(PETSC_USE_COMPLEX)
582:   if (ki) SETERRQ(PETSC_COMM_SELF,1,"Not implemented for complex eigenvalues with real scalars");
583: #endif
584:   NEPComputeResidualNorm_Private(nep,kr,xr,norm);
585:   VecDestroy(&xr);
586:   VecDestroy(&xi);
587:   return(0);
588: }

592: /*
593:    NEPComputeRelativeError_Private - Computes the relative error bound
594:    associated with an eigenpair.
595: */
596: PetscErrorCode NEPComputeRelativeError_Private(NEP nep,PetscScalar lambda,Vec x,PetscReal *error)
597: {
599:   PetscReal      norm,er;

602:   NEPComputeResidualNorm_Private(nep,lambda,x,&norm);
603:   VecNorm(x,NORM_2,&er);
604:   if (PetscAbsScalar(lambda) > norm) {
605:     *error = norm/(PetscAbsScalar(lambda)*er);
606:   } else {
607:     *error = norm/er;
608:   }
609:   return(0);
610: }

614: /*@
615:    NEPComputeRelativeError - Computes the relative error bound associated
616:    with the i-th computed eigenpair.

618:    Collective on NEP

620:    Input Parameter:
621: +  nep - the nonlinear eigensolver context
622: -  i   - the solution index

624:    Output Parameter:
625: .  error - the relative error bound, computed as ||T(lambda)x||_2/||lambda*x||_2
626:    where lambda is the eigenvalue and x is the eigenvector.
627:    If lambda=0 the relative error is computed as ||T(lambda)x||_2/||x||_2.

629:    Level: beginner

631: .seealso: NEPSolve(), NEPComputeResidualNorm(), NEPGetErrorEstimate()
632: @*/
633: PetscErrorCode NEPComputeRelativeError(NEP nep,PetscInt i,PetscReal *error)
634: {
636:   Vec            xr,xi;
637:   PetscScalar    kr,ki;

643:   NEPCheckSolved(nep,1);
644:   BVGetVec(nep->V,&xr);
645:   BVGetVec(nep->V,&xi);
646:   NEPGetEigenpair(nep,i,&kr,&ki,xr,xi);
647: #if !defined(PETSC_USE_COMPLEX)
648:   if (ki) SETERRQ(PETSC_COMM_SELF,1,"Not implemented for complex eigenvalues with real scalars");
649: #endif
650:   NEPComputeRelativeError_Private(nep,kr,xr,error);
651:   VecDestroy(&xr);
652:   VecDestroy(&xi);
653:   return(0);
654: }

658: /*@
659:    NEPComputeFunction - Computes the function matrix T(lambda) that has been
660:    set with NEPSetFunction().

662:    Collective on NEP and Mat

664:    Input Parameters:
665: +  nep    - the NEP context
666: -  lambda - the scalar argument

668:    Output Parameters:
669: +  A   - Function matrix
670: -  B   - optional preconditioning matrix

672:    Notes:
673:    NEPComputeFunction() is typically used within nonlinear eigensolvers
674:    implementations, so most users would not generally call this routine
675:    themselves.

677:    Level: developer

679: .seealso: NEPSetFunction(), NEPGetFunction()
680: @*/
681: PetscErrorCode NEPComputeFunction(NEP nep,PetscScalar lambda,Mat A,Mat B)
682: {
684:   PetscInt       i;
685:   PetscScalar    alpha;


690:   if (nep->split) {

692:     MatZeroEntries(A);
693:     for (i=0;i<nep->nt;i++) {
694:       FNEvaluateFunction(nep->f[i],lambda,&alpha);
695:       MatAXPY(A,alpha,nep->A[i],nep->mstr);
696:     }
697:     if (A != B) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Not implemented");

699:   } else {

701:     if (!nep->computefunction) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_USER,"Must call NEPSetFunction() first");

703:     PetscLogEventBegin(NEP_FunctionEval,nep,A,B,0);

705:     PetscStackPush("NEP user Function function");
706:     (*nep->computefunction)(nep,lambda,A,B,nep->functionctx);
707:     PetscStackPop;

709:     PetscLogEventEnd(NEP_FunctionEval,nep,A,B,0);
710:     nep->nfuncs++;

712:   }
713:   return(0);
714: }

718: /*@
719:    NEPComputeJacobian - Computes the Jacobian matrix T'(lambda) that has been
720:    set with NEPSetJacobian().

722:    Collective on NEP and Mat

724:    Input Parameters:
725: +  nep    - the NEP context
726: -  lambda - the scalar argument

728:    Output Parameters:
729: .  A   - Jacobian matrix

731:    Notes:
732:    Most users should not need to explicitly call this routine, as it
733:    is used internally within the nonlinear eigensolvers.

735:    Level: developer

737: .seealso: NEPSetJacobian(), NEPGetJacobian()
738: @*/
739: PetscErrorCode NEPComputeJacobian(NEP nep,PetscScalar lambda,Mat A)
740: {
742:   PetscInt       i;
743:   PetscScalar    alpha;


748:   if (nep->split) {

750:     MatZeroEntries(A);
751:     for (i=0;i<nep->nt;i++) {
752:       FNEvaluateDerivative(nep->f[i],lambda,&alpha);
753:       MatAXPY(A,alpha,nep->A[i],nep->mstr);
754:     }

756:   } else {

758:     if (!nep->computejacobian) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_USER,"Must call NEPSetJacobian() first");

760:     PetscLogEventBegin(NEP_JacobianEval,nep,A,0,0);

762:     PetscStackPush("NEP user Jacobian function");
763:     (*nep->computejacobian)(nep,lambda,A,nep->jacobianctx);
764:     PetscStackPop;

766:     PetscLogEventEnd(NEP_JacobianEval,nep,A,0,0);

768:   }
769:   return(0);
770: }