Actual source code: nepsolve.c

slepc-3.15.0 2021-03-31
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:   "   note = \"To appear\",\n"
 31:   "   url = \"https://arxiv.org/abs/1910.11712\",\n"
 32:   "   year = \"2021\"\n"
 33:   "}\n";

 35: PetscErrorCode NEPComputeVectors(NEP nep)
 36: {

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

 48: /*@
 49:    NEPSolve - Solves the nonlinear eigensystem.

 51:    Collective on nep

 53:    Input Parameter:
 54: .  nep - eigensolver context obtained from NEPCreate()

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

 64:    Level: beginner

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

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

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

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

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

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

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

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

118:   /* Remove the initial subspace */
119:   nep->nini = 0;

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

126: /*@
127:    NEPProjectOperator - Computes the projection of the nonlinear operator.

129:    Collective on nep

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

136:    Notes:
137:    This is available for split operator only.

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

146:    Level: developer

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

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

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

174:    Collective on nep

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

182:    Output Parameters:
183: +  y   - result vector
184: .  A   - Function matrix
185: -  B   - optional preconditioning matrix

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

194:    Level: developer

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


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

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

230:    Collective on nep

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

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

243:    Level: developer

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


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

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

282:    Collective on nep

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

290:    Output Parameters:
291: +  y   - result vector
292: -  A   - Jacobian matrix

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

301:    Level: developer

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


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

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

338:    Not Collective

340:    Input Parameter:
341: .  nep - the nonlinear eigensolver context

343:    Output Parameter:
344: .  its - number of iterations

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

353:    Level: intermediate

355: .seealso: NEPGetConvergedReason(), NEPSetTolerances()
356: @*/
357: PetscErrorCode NEPGetIterationNumber(NEP nep,PetscInt *its)
358: {
362:   *its = nep->its;
363:   return(0);
364: }

366: /*@
367:    NEPGetConverged - Gets the number of converged eigenpairs.

369:    Not Collective

371:    Input Parameter:
372: .  nep - the nonlinear eigensolver context

374:    Output Parameter:
375: .  nconv - number of converged eigenpairs

377:    Note:
378:    This function should be called after NEPSolve() has finished.

380:    Level: beginner

382: .seealso: NEPSetDimensions(), NEPSolve()
383: @*/
384: PetscErrorCode NEPGetConverged(NEP nep,PetscInt *nconv)
385: {
389:   NEPCheckSolved(nep,1);
390:   *nconv = nep->nconv;
391:   return(0);
392: }

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

398:    Not Collective

400:    Input Parameter:
401: .  nep - the nonlinear eigensolver context

403:    Output Parameter:
404: .  reason - negative value indicates diverged, positive value converged

406:    Notes:

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

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

419:    Level: intermediate

421: .seealso: NEPSetTolerances(), NEPSolve(), NEPConvergedReason
422: @*/
423: PetscErrorCode NEPGetConvergedReason(NEP nep,NEPConvergedReason *reason)
424: {
428:   NEPCheckSolved(nep,1);
429:   *reason = nep->reason;
430:   return(0);
431: }

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

437:    Logically Collective on nep

439:    Input Parameters:
440: +  nep - nonlinear eigensolver context
441: -  i   - index of the solution

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

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

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

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

464:    Level: beginner

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

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

482:   NEPComputeVectors(nep);
483:   k = nep->perm[i];

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

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

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

502:    Logically Collective on nep

504:    Input Parameters:
505: +  nep - eigensolver context
506: -  i   - index of the solution

508:    Output Parameters:
509: +  Wr   - real part of left eigenvector
510: -  Wi   - imaginary part of left eigenvector

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

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

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

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

528:    Level: intermediate

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

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

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

555:    Not Collective

557:    Input Parameter:
558: +  nep - nonlinear eigensolver context
559: -  i   - index of eigenpair

561:    Output Parameter:
562: .  errest - the error estimate

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

568:    Level: advanced

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

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

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

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

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

614:    Collective on nep

616:    Input Parameter:
617: +  nep  - the nonlinear eigensolver context
618: .  i    - the solution index
619: -  type - the type of error to compute

621:    Output Parameter:
622: .  error - the error

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

629:    Level: beginner

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

647:   NEPCheckSolved(nep,1);

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

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

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

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

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

720:    Collective on nep

722:    Input Parameters:
723: +  nep    - the NEP context
724: -  lambda - the scalar argument

726:    Output Parameters:
727: +  A   - Function matrix
728: -  B   - optional preconditioning matrix

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

735:    Level: developer

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

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

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

773:    Collective on nep

775:    Input Parameters:
776: +  nep    - the NEP context
777: -  lambda - the scalar argument

779:    Output Parameters:
780: .  A   - Jacobian matrix

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

786:    Level: developer

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

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