Actual source code: nepsolve.c

slepc-3.16.1 2021-11-17
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   - Function matrix
187: -  B   - optional 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:     NEPComputeFunction(nep,lambda,A,B);
224:     MatMult(A,x,y);
225:   }
226:   return(0);
227: }

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

232:    Collective on nep

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

240:    Output Parameters:
241: +  y   - result vector
242: .  A   - Function matrix
243: -  B   - optional preconditioning matrix

245:    Level: developer

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


264:   VecConjugate(x);
265:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
266:     VecSet(y,0.0);
267:     for (i=0;i<nep->nt;i++) {
268:       FNEvaluateFunction(nep->f[i],lambda,&alpha);
269:       MatMultTranspose(nep->A[i],x,v);
270:       VecAXPY(y,alpha,v);
271:     }
272:   } else {
273:     NEPComputeFunction(nep,lambda,A,B);
274:     MatMultTranspose(A,x,y);
275:   }
276:   VecConjugate(x);
277:   VecConjugate(y);
278:   return(0);
279: }

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

284:    Collective on nep

286:    Input Parameters:
287: +  nep    - the nonlinear eigensolver context
288: .  lambda - scalar argument
289: .  x      - vector to be multiplied against
290: -  v      - workspace vector (used only in the case of split form)

292:    Output Parameters:
293: +  y   - result vector
294: -  A   - Jacobian matrix

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

303:    Level: developer

305: .seealso: NEPSetSplitOperator(), NEPComputeJacobian()
306: @*/
307: PetscErrorCode NEPApplyJacobian(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A)
308: {
310:   PetscInt       i;
311:   PetscScalar    alpha;


321:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
322:     VecSet(y,0.0);
323:     for (i=0;i<nep->nt;i++) {
324:       FNEvaluateDerivative(nep->f[i],lambda,&alpha);
325:       MatMult(nep->A[i],x,v);
326:       VecAXPY(y,alpha,v);
327:     }
328:   } else {
329:     NEPComputeJacobian(nep,lambda,A);
330:     MatMult(A,x,y);
331:   }
332:   return(0);
333: }

335: /*@
336:    NEPGetIterationNumber - Gets the current iteration number. If the
337:    call to NEPSolve() is complete, then it returns the number of iterations
338:    carried out by the solution method.

340:    Not Collective

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

345:    Output Parameter:
346: .  its - number of iterations

348:    Note:
349:    During the i-th iteration this call returns i-1. If NEPSolve() is
350:    complete, then parameter "its" contains either the iteration number at
351:    which convergence was successfully reached, or failure was detected.
352:    Call NEPGetConvergedReason() to determine if the solver converged or
353:    failed and why.

355:    Level: intermediate

357: .seealso: NEPGetConvergedReason(), NEPSetTolerances()
358: @*/
359: PetscErrorCode NEPGetIterationNumber(NEP nep,PetscInt *its)
360: {
364:   *its = nep->its;
365:   return(0);
366: }

368: /*@
369:    NEPGetConverged - Gets the number of converged eigenpairs.

371:    Not Collective

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

376:    Output Parameter:
377: .  nconv - number of converged eigenpairs

379:    Note:
380:    This function should be called after NEPSolve() has finished.

382:    Level: beginner

384: .seealso: NEPSetDimensions(), NEPSolve()
385: @*/
386: PetscErrorCode NEPGetConverged(NEP nep,PetscInt *nconv)
387: {
391:   NEPCheckSolved(nep,1);
392:   *nconv = nep->nconv;
393:   return(0);
394: }

396: /*@
397:    NEPGetConvergedReason - Gets the reason why the NEPSolve() iteration was
398:    stopped.

400:    Not Collective

402:    Input Parameter:
403: .  nep - the nonlinear eigensolver context

405:    Output Parameter:
406: .  reason - negative value indicates diverged, positive value converged

408:    Notes:

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

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

421:    Level: intermediate

423: .seealso: NEPSetTolerances(), NEPSolve(), NEPConvergedReason
424: @*/
425: PetscErrorCode NEPGetConvergedReason(NEP nep,NEPConvergedReason *reason)
426: {
430:   NEPCheckSolved(nep,1);
431:   *reason = nep->reason;
432:   return(0);
433: }

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

439:    Logically Collective on nep

441:    Input Parameters:
442: +  nep - nonlinear eigensolver context
443: -  i   - index of the solution

445:    Output Parameters:
446: +  eigr - real part of eigenvalue
447: .  eigi - imaginary part of eigenvalue
448: .  Vr   - real part of eigenvector
449: -  Vi   - imaginary part of eigenvector

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

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

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

466:    Level: beginner

468: .seealso: NEPSolve(), NEPGetConverged(), NEPSetWhichEigenpairs(), NEPGetLeftEigenvector()
469: @*/
470: PetscErrorCode NEPGetEigenpair(NEP nep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
471: {
472:   PetscInt       k;

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

484:   NEPComputeVectors(nep);
485:   k = nep->perm[i];

487:   /* eigenvalue */
488: #if defined(PETSC_USE_COMPLEX)
489:   if (eigr) *eigr = nep->eigr[k];
490:   if (eigi) *eigi = 0;
491: #else
492:   if (eigr) *eigr = nep->eigr[k];
493:   if (eigi) *eigi = nep->eigi[k];
494: #endif

496:   /* eigenvector */
497:   BV_GetEigenvector(nep->V,k,nep->eigi[k],Vr,Vi);
498:   return(0);
499: }

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

504:    Logically Collective on nep

506:    Input Parameters:
507: +  nep - eigensolver context
508: -  i   - index of the solution

510:    Output Parameters:
511: +  Wr   - real part of left eigenvector
512: -  Wi   - imaginary part of left eigenvector

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

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

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

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

530:    Level: intermediate

532: .seealso: NEPGetEigenpair(), NEPGetConverged(), NEPSetWhichEigenpairs(), NEPSetTwoSided()
533: @*/
534: PetscErrorCode NEPGetLeftEigenvector(NEP nep,PetscInt i,Vec Wr,Vec Wi)
535: {
537:   PetscInt       k;

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

553: /*@
554:    NEPGetErrorEstimate - Returns the error estimate associated to the i-th
555:    computed eigenpair.

557:    Not Collective

559:    Input Parameters:
560: +  nep - nonlinear eigensolver context
561: -  i   - index of eigenpair

563:    Output Parameter:
564: .  errest - the error estimate

566:    Notes:
567:    This is the error estimate used internally by the eigensolver. The actual
568:    error bound can be computed with NEPComputeError().

570:    Level: advanced

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

585: /*
586:    NEPComputeResidualNorm_Private - Computes the norm of the residual vector
587:    associated with an eigenpair.

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

601:   y = w[0];
602:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) z = w[1];
603:   if (adj) {
604:     NEPApplyAdjoint(nep,lambda,x,z,y,nep->function,nep->function_pre);
605:   } else {
606:     NEPApplyFunction(nep,lambda,x,z,y,nep->function,nep->function_pre);
607:   }
608:   VecNorm(y,NORM_2,norm);
609:   return(0);
610: }

612: /*@
613:    NEPComputeError - Computes the error (based on the residual norm) associated
614:    with the i-th computed eigenpair.

616:    Collective on nep

618:    Input Parameters:
619: +  nep  - the nonlinear eigensolver context
620: .  i    - the solution index
621: -  type - the type of error to compute

623:    Output Parameter:
624: .  error - the error

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

631:    Level: beginner

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

649:   NEPCheckSolved(nep,1);

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

667:   /* compute residual norms */
668:   NEPGetEigenpair(nep,i,&kr,&ki,xr,xi);
669: #if !defined(PETSC_USE_COMPLEX)
670:   if (ki) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Not implemented for complex eigenvalues with real scalars");
671: #endif
672:   NEPComputeResidualNorm_Private(nep,PETSC_FALSE,kr,xr,nep->work,error);
673:   VecNorm(xr,NORM_2,&er);

675:   /* if two-sided, compute left residual norm and take the maximum */
676:   if (nep->twosided) {
677:     NEPGetLeftEigenvector(nep,i,xr,xi);
678:     NEPComputeResidualNorm_Private(nep,PETSC_TRUE,kr,xr,nep->work,&errorl);
679:     *error = PetscMax(*error,errorl);
680:   }

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

718: /*@
719:    NEPComputeFunction - Computes the function matrix T(lambda) that has been
720:    set with NEPSetFunction().

722:    Collective on nep

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

728:    Output Parameters:
729: +  A   - Function matrix
730: -  B   - optional preconditioning matrix

732:    Notes:
733:    NEPComputeFunction() is typically used within nonlinear eigensolvers
734:    implementations, so most users would not generally call this routine
735:    themselves.

737:    Level: developer

739: .seealso: NEPSetFunction(), NEPGetFunction()
740: @*/
741: PetscErrorCode NEPComputeFunction(NEP nep,PetscScalar lambda,Mat A,Mat B)
742: {
744:   PetscInt       i;
745:   PetscScalar    alpha;

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

771: /*@
772:    NEPComputeJacobian - Computes the Jacobian matrix T'(lambda) that has been
773:    set with NEPSetJacobian().

775:    Collective on nep

777:    Input Parameters:
778: +  nep    - the NEP context
779: -  lambda - the scalar argument

781:    Output Parameters:
782: .  A   - Jacobian matrix

784:    Notes:
785:    Most users should not need to explicitly call this routine, as it
786:    is used internally within the nonlinear eigensolvers.

788:    Level: developer

790: .seealso: NEPSetJacobian(), NEPGetJacobian()
791: @*/
792: PetscErrorCode NEPComputeJacobian(NEP nep,PetscScalar lambda,Mat A)
793: {
795:   PetscInt       i;
796:   PetscScalar    alpha;

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