Actual source code: nepbasic.c

slepc-3.15.2 2021-09-20
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:    Basic NEP routines
 12: */

 14: #include <slepc/private/nepimpl.h>

 16: /* Logging support */
 17: PetscClassId      NEP_CLASSID = 0;
 18: PetscLogEvent     NEP_SetUp = 0,NEP_Solve = 0,NEP_Refine = 0,NEP_FunctionEval = 0,NEP_JacobianEval = 0,NEP_Resolvent = 0;

 20: /* List of registered NEP routines */
 21: PetscFunctionList NEPList = 0;
 22: PetscBool         NEPRegisterAllCalled = PETSC_FALSE;

 24: /* List of registered NEP monitors */
 25: PetscFunctionList NEPMonitorList              = NULL;
 26: PetscFunctionList NEPMonitorCreateList        = NULL;
 27: PetscFunctionList NEPMonitorDestroyList       = NULL;
 28: PetscBool         NEPMonitorRegisterAllCalled = PETSC_FALSE;

 30: /*@
 31:    NEPCreate - Creates the default NEP context.

 33:    Collective

 35:    Input Parameter:
 36: .  comm - MPI communicator

 38:    Output Parameter:
 39: .  nep - location to put the NEP context

 41:    Level: beginner

 43: .seealso: NEPSetUp(), NEPSolve(), NEPDestroy(), NEP
 44: @*/
 45: PetscErrorCode NEPCreate(MPI_Comm comm,NEP *outnep)
 46: {
 48:   NEP            nep;

 52:   *outnep = 0;
 53:   NEPInitializePackage();
 54:   SlepcHeaderCreate(nep,NEP_CLASSID,"NEP","Nonlinear Eigenvalue Problem","NEP",comm,NEPDestroy,NEPView);

 56:   nep->max_it          = PETSC_DEFAULT;
 57:   nep->nev             = 1;
 58:   nep->ncv             = PETSC_DEFAULT;
 59:   nep->mpd             = PETSC_DEFAULT;
 60:   nep->nini            = 0;
 61:   nep->target          = 0.0;
 62:   nep->tol             = PETSC_DEFAULT;
 63:   nep->conv            = NEP_CONV_NORM;
 64:   nep->stop            = NEP_STOP_BASIC;
 65:   nep->which           = (NEPWhich)0;
 66:   nep->problem_type    = (NEPProblemType)0;
 67:   nep->refine          = NEP_REFINE_NONE;
 68:   nep->npart           = 1;
 69:   nep->rtol            = PETSC_DEFAULT;
 70:   nep->rits            = PETSC_DEFAULT;
 71:   nep->scheme          = (NEPRefineScheme)0;
 72:   nep->trackall        = PETSC_FALSE;
 73:   nep->twosided        = PETSC_FALSE;

 75:   nep->computefunction = NULL;
 76:   nep->computejacobian = NULL;
 77:   nep->functionctx     = NULL;
 78:   nep->jacobianctx     = NULL;
 79:   nep->converged       = NEPConvergedRelative;
 80:   nep->convergeduser   = NULL;
 81:   nep->convergeddestroy= NULL;
 82:   nep->stopping        = NEPStoppingBasic;
 83:   nep->stoppinguser    = NULL;
 84:   nep->stoppingdestroy = NULL;
 85:   nep->convergedctx    = NULL;
 86:   nep->stoppingctx     = NULL;
 87:   nep->numbermonitors  = 0;

 89:   nep->ds              = NULL;
 90:   nep->V               = NULL;
 91:   nep->W               = NULL;
 92:   nep->rg              = NULL;
 93:   nep->function        = NULL;
 94:   nep->function_pre    = NULL;
 95:   nep->jacobian        = NULL;
 96:   nep->A               = NULL;
 97:   nep->f               = NULL;
 98:   nep->nt              = 0;
 99:   nep->mstr            = UNKNOWN_NONZERO_PATTERN;
100:   nep->IS              = NULL;
101:   nep->eigr            = NULL;
102:   nep->eigi            = NULL;
103:   nep->errest          = NULL;
104:   nep->perm            = NULL;
105:   nep->nwork           = 0;
106:   nep->work            = NULL;
107:   nep->data            = NULL;

109:   nep->state           = NEP_STATE_INITIAL;
110:   nep->nconv           = 0;
111:   nep->its             = 0;
112:   nep->n               = 0;
113:   nep->nloc            = 0;
114:   nep->nrma            = NULL;
115:   nep->fui             = (NEPUserInterface)0;
116:   nep->useds           = PETSC_FALSE;
117:   nep->resolvent       = NULL;
118:   nep->reason          = NEP_CONVERGED_ITERATING;

120:   PetscNewLog(nep,&nep->sc);
121:   *outnep = nep;
122:   return(0);
123: }

125: /*@C
126:    NEPSetType - Selects the particular solver to be used in the NEP object.

128:    Logically Collective on nep

130:    Input Parameters:
131: +  nep      - the nonlinear eigensolver context
132: -  type     - a known method

134:    Options Database Key:
135: .  -nep_type <method> - Sets the method; use -help for a list
136:     of available methods

138:    Notes:
139:    See "slepc/include/slepcnep.h" for available methods.

141:    Normally, it is best to use the NEPSetFromOptions() command and
142:    then set the NEP type from the options database rather than by using
143:    this routine.  Using the options database provides the user with
144:    maximum flexibility in evaluating the different available methods.
145:    The NEPSetType() routine is provided for those situations where it
146:    is necessary to set the iterative solver independently of the command
147:    line or options database.

149:    Level: intermediate

151: .seealso: NEPType
152: @*/
153: PetscErrorCode NEPSetType(NEP nep,NEPType type)
154: {
155:   PetscErrorCode ierr,(*r)(NEP);
156:   PetscBool      match;


162:   PetscObjectTypeCompare((PetscObject)nep,type,&match);
163:   if (match) return(0);

165:   PetscFunctionListFind(NEPList,type,&r);
166:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown NEP type given: %s",type);

168:   if (nep->ops->destroy) { (*nep->ops->destroy)(nep); }
169:   PetscMemzero(nep->ops,sizeof(struct _NEPOps));

171:   nep->state = NEP_STATE_INITIAL;
172:   PetscObjectChangeTypeName((PetscObject)nep,type);
173:   (*r)(nep);
174:   return(0);
175: }

177: /*@C
178:    NEPGetType - Gets the NEP type as a string from the NEP object.

180:    Not Collective

182:    Input Parameter:
183: .  nep - the eigensolver context

185:    Output Parameter:
186: .  name - name of NEP method

188:    Level: intermediate

190: .seealso: NEPSetType()
191: @*/
192: PetscErrorCode NEPGetType(NEP nep,NEPType *type)
193: {
197:   *type = ((PetscObject)nep)->type_name;
198:   return(0);
199: }

201: /*@C
202:    NEPRegister - Adds a method to the nonlinear eigenproblem solver package.

204:    Not Collective

206:    Input Parameters:
207: +  name - name of a new user-defined solver
208: -  function - routine to create the solver context

210:    Notes:
211:    NEPRegister() may be called multiple times to add several user-defined solvers.

213:    Sample usage:
214: .vb
215:     NEPRegister("my_solver",MySolverCreate);
216: .ve

218:    Then, your solver can be chosen with the procedural interface via
219: $     NEPSetType(nep,"my_solver")
220:    or at runtime via the option
221: $     -nep_type my_solver

223:    Level: advanced

225: .seealso: NEPRegisterAll()
226: @*/
227: PetscErrorCode NEPRegister(const char *name,PetscErrorCode (*function)(NEP))
228: {

232:   NEPInitializePackage();
233:   PetscFunctionListAdd(&NEPList,name,function);
234:   return(0);
235: }

237: /*@C
238:    NEPMonitorRegister - Adds NEP monitor routine.

240:    Not Collective

242:    Input Parameters:
243: +  name    - name of a new monitor routine
244: .  vtype   - a PetscViewerType for the output
245: .  format  - a PetscViewerFormat for the output
246: .  monitor - monitor routine
247: .  create  - creation routine, or NULL
248: -  destroy - destruction routine, or NULL

250:    Notes:
251:    NEPMonitorRegister() may be called multiple times to add several user-defined monitors.

253:    Sample usage:
254: .vb
255:    NEPMonitorRegister("my_monitor",PETSCVIEWERASCII,PETSC_VIEWER_ASCII_INFO_DETAIL,MyMonitor,NULL,NULL);
256: .ve

258:    Then, your monitor can be chosen with the procedural interface via
259: $      NEPMonitorSetFromOptions(nep,"-nep_monitor_my_monitor","my_monitor",NULL)
260:    or at runtime via the option
261: $      -nep_monitor_my_monitor

263:    Level: advanced

265: .seealso: NEPMonitorRegisterAll()
266: @*/
267: PetscErrorCode NEPMonitorRegister(const char name[],PetscViewerType vtype,PetscViewerFormat format,PetscErrorCode (*monitor)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode (*create)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode (*destroy)(PetscViewerAndFormat**))
268: {
269:   char           key[PETSC_MAX_PATH_LEN];

273:   NEPInitializePackage();
274:   SlepcMonitorMakeKey_Internal(name,vtype,format,key);
275:   PetscFunctionListAdd(&NEPMonitorList,key,monitor);
276:   if (create)  { PetscFunctionListAdd(&NEPMonitorCreateList,key,create); }
277:   if (destroy) { PetscFunctionListAdd(&NEPMonitorDestroyList,key,destroy); }
278:   return(0);
279: }

281: /*
282:    NEPReset_Problem - Destroys the problem matrices.
283: @*/
284: PetscErrorCode NEPReset_Problem(NEP nep)
285: {
287:   PetscInt       i;

291:   MatDestroy(&nep->function);
292:   MatDestroy(&nep->function_pre);
293:   MatDestroy(&nep->jacobian);
294:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
295:     MatDestroyMatrices(nep->nt,&nep->A);
296:     for (i=0;i<nep->nt;i++) {
297:       FNDestroy(&nep->f[i]);
298:     }
299:     PetscFree(nep->f);
300:     PetscFree(nep->nrma);
301:     nep->nt = 0;
302:   }
303:   return(0);
304: }
305: /*@
306:    NEPReset - Resets the NEP context to the initial state (prior to setup)
307:    and destroys any allocated Vecs and Mats.

309:    Collective on nep

311:    Input Parameter:
312: .  nep - eigensolver context obtained from NEPCreate()

314:    Level: advanced

316: .seealso: NEPDestroy()
317: @*/
318: PetscErrorCode NEPReset(NEP nep)
319: {

324:   if (!nep) return(0);
325:   if (nep->ops->reset) { (nep->ops->reset)(nep); }
326:   if (nep->refineksp) { KSPReset(nep->refineksp); }
327:   NEPReset_Problem(nep);
328:   BVDestroy(&nep->V);
329:   BVDestroy(&nep->W);
330:   VecDestroyVecs(nep->nwork,&nep->work);
331:   MatDestroy(&nep->resolvent);
332:   nep->nwork = 0;
333:   nep->state = NEP_STATE_INITIAL;
334:   return(0);
335: }

337: /*@C
338:    NEPDestroy - Destroys the NEP context.

340:    Collective on nep

342:    Input Parameter:
343: .  nep - eigensolver context obtained from NEPCreate()

345:    Level: beginner

347: .seealso: NEPCreate(), NEPSetUp(), NEPSolve()
348: @*/
349: PetscErrorCode NEPDestroy(NEP *nep)
350: {

354:   if (!*nep) return(0);
356:   if (--((PetscObject)(*nep))->refct > 0) { *nep = 0; return(0); }
357:   NEPReset(*nep);
358:   if ((*nep)->ops->destroy) { (*(*nep)->ops->destroy)(*nep); }
359:   if ((*nep)->eigr) {
360:     PetscFree4((*nep)->eigr,(*nep)->eigi,(*nep)->errest,(*nep)->perm);
361:   }
362:   RGDestroy(&(*nep)->rg);
363:   DSDestroy(&(*nep)->ds);
364:   KSPDestroy(&(*nep)->refineksp);
365:   PetscSubcommDestroy(&(*nep)->refinesubc);
366:   PetscFree((*nep)->sc);
367:   /* just in case the initial vectors have not been used */
368:   SlepcBasisDestroy_Private(&(*nep)->nini,&(*nep)->IS);
369:   if ((*nep)->convergeddestroy) {
370:     (*(*nep)->convergeddestroy)((*nep)->convergedctx);
371:   }
372:   NEPMonitorCancel(*nep);
373:   PetscHeaderDestroy(nep);
374:   return(0);
375: }

377: /*@
378:    NEPSetBV - Associates a basis vectors object to the nonlinear eigensolver.

380:    Collective on nep

382:    Input Parameters:
383: +  nep - eigensolver context obtained from NEPCreate()
384: -  bv  - the basis vectors object

386:    Note:
387:    Use NEPGetBV() to retrieve the basis vectors context (for example,
388:    to free it at the end of the computations).

390:    Level: advanced

392: .seealso: NEPGetBV()
393: @*/
394: PetscErrorCode NEPSetBV(NEP nep,BV bv)
395: {

402:   PetscObjectReference((PetscObject)bv);
403:   BVDestroy(&nep->V);
404:   nep->V = bv;
405:   PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->V);
406:   return(0);
407: }

409: /*@
410:    NEPGetBV - Obtain the basis vectors object associated to the nonlinear
411:    eigensolver object.

413:    Not Collective

415:    Input Parameters:
416: .  nep - eigensolver context obtained from NEPCreate()

418:    Output Parameter:
419: .  bv - basis vectors context

421:    Level: advanced

423: .seealso: NEPSetBV()
424: @*/
425: PetscErrorCode NEPGetBV(NEP nep,BV *bv)
426: {

432:   if (!nep->V) {
433:     BVCreate(PetscObjectComm((PetscObject)nep),&nep->V);
434:     PetscObjectIncrementTabLevel((PetscObject)nep->V,(PetscObject)nep,0);
435:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->V);
436:     PetscObjectSetOptions((PetscObject)nep->V,((PetscObject)nep)->options);
437:   }
438:   *bv = nep->V;
439:   return(0);
440: }

442: /*@
443:    NEPSetRG - Associates a region object to the nonlinear eigensolver.

445:    Collective on nep

447:    Input Parameters:
448: +  nep - eigensolver context obtained from NEPCreate()
449: -  rg  - the region object

451:    Note:
452:    Use NEPGetRG() to retrieve the region context (for example,
453:    to free it at the end of the computations).

455:    Level: advanced

457: .seealso: NEPGetRG()
458: @*/
459: PetscErrorCode NEPSetRG(NEP nep,RG rg)
460: {

467:   PetscObjectReference((PetscObject)rg);
468:   RGDestroy(&nep->rg);
469:   nep->rg = rg;
470:   PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->rg);
471:   return(0);
472: }

474: /*@
475:    NEPGetRG - Obtain the region object associated to the
476:    nonlinear eigensolver object.

478:    Not Collective

480:    Input Parameters:
481: .  nep - eigensolver context obtained from NEPCreate()

483:    Output Parameter:
484: .  rg - region context

486:    Level: advanced

488: .seealso: NEPSetRG()
489: @*/
490: PetscErrorCode NEPGetRG(NEP nep,RG *rg)
491: {

497:   if (!nep->rg) {
498:     RGCreate(PetscObjectComm((PetscObject)nep),&nep->rg);
499:     PetscObjectIncrementTabLevel((PetscObject)nep->rg,(PetscObject)nep,0);
500:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->rg);
501:     PetscObjectSetOptions((PetscObject)nep->rg,((PetscObject)nep)->options);
502:   }
503:   *rg = nep->rg;
504:   return(0);
505: }

507: /*@
508:    NEPSetDS - Associates a direct solver object to the nonlinear eigensolver.

510:    Collective on nep

512:    Input Parameters:
513: +  nep - eigensolver context obtained from NEPCreate()
514: -  ds  - the direct solver object

516:    Note:
517:    Use NEPGetDS() to retrieve the direct solver context (for example,
518:    to free it at the end of the computations).

520:    Level: advanced

522: .seealso: NEPGetDS()
523: @*/
524: PetscErrorCode NEPSetDS(NEP nep,DS ds)
525: {

532:   PetscObjectReference((PetscObject)ds);
533:   DSDestroy(&nep->ds);
534:   nep->ds = ds;
535:   PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->ds);
536:   return(0);
537: }

539: /*@
540:    NEPGetDS - Obtain the direct solver object associated to the
541:    nonlinear eigensolver object.

543:    Not Collective

545:    Input Parameters:
546: .  nep - eigensolver context obtained from NEPCreate()

548:    Output Parameter:
549: .  ds - direct solver context

551:    Level: advanced

553: .seealso: NEPSetDS()
554: @*/
555: PetscErrorCode NEPGetDS(NEP nep,DS *ds)
556: {

562:   if (!nep->ds) {
563:     DSCreate(PetscObjectComm((PetscObject)nep),&nep->ds);
564:     PetscObjectIncrementTabLevel((PetscObject)nep->ds,(PetscObject)nep,0);
565:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->ds);
566:     PetscObjectSetOptions((PetscObject)nep->ds,((PetscObject)nep)->options);
567:   }
568:   *ds = nep->ds;
569:   return(0);
570: }

572: /*@
573:    NEPRefineGetKSP - Obtain the ksp object used by the eigensolver
574:    object in the refinement phase.

576:    Not Collective

578:    Input Parameters:
579: .  nep - eigensolver context obtained from NEPCreate()

581:    Output Parameter:
582: .  ksp - ksp context

584:    Level: advanced

586: .seealso: NEPSetRefine()
587: @*/
588: PetscErrorCode NEPRefineGetKSP(NEP nep,KSP *ksp)
589: {

595:   if (!nep->refineksp) {
596:     if (nep->npart>1) {
597:       /* Split in subcomunicators */
598:       PetscSubcommCreate(PetscObjectComm((PetscObject)nep),&nep->refinesubc);
599:       PetscSubcommSetNumber(nep->refinesubc,nep->npart);
600:       PetscSubcommSetType(nep->refinesubc,PETSC_SUBCOMM_CONTIGUOUS);
601:       PetscLogObjectMemory((PetscObject)nep,sizeof(PetscSubcomm));
602:     }
603:     KSPCreate((nep->npart==1)?PetscObjectComm((PetscObject)nep):PetscSubcommChild(nep->refinesubc),&nep->refineksp);
604:     PetscObjectIncrementTabLevel((PetscObject)nep->refineksp,(PetscObject)nep,0);
605:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->refineksp);
606:     PetscObjectSetOptions((PetscObject)nep->refineksp,((PetscObject)nep)->options);
607:     KSPSetOptionsPrefix(*ksp,((PetscObject)nep)->prefix);
608:     KSPAppendOptionsPrefix(*ksp,"nep_refine_");
609:     KSPSetTolerances(nep->refineksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
610:   }
611:   *ksp = nep->refineksp;
612:   return(0);
613: }

615: /*@
616:    NEPSetTarget - Sets the value of the target.

618:    Logically Collective on nep

620:    Input Parameters:
621: +  nep    - eigensolver context
622: -  target - the value of the target

624:    Options Database Key:
625: .  -nep_target <scalar> - the value of the target

627:    Notes:
628:    The target is a scalar value used to determine the portion of the spectrum
629:    of interest. It is used in combination with NEPSetWhichEigenpairs().

631:    In the case of complex scalars, a complex value can be provided in the
632:    command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
633:    -nep_target 1.0+2.0i

635:    Level: intermediate

637: .seealso: NEPGetTarget(), NEPSetWhichEigenpairs()
638: @*/
639: PetscErrorCode NEPSetTarget(NEP nep,PetscScalar target)
640: {
644:   nep->target = target;
645:   return(0);
646: }

648: /*@
649:    NEPGetTarget - Gets the value of the target.

651:    Not Collective

653:    Input Parameter:
654: .  nep - eigensolver context

656:    Output Parameter:
657: .  target - the value of the target

659:    Note:
660:    If the target was not set by the user, then zero is returned.

662:    Level: intermediate

664: .seealso: NEPSetTarget()
665: @*/
666: PetscErrorCode NEPGetTarget(NEP nep,PetscScalar* target)
667: {
671:   *target = nep->target;
672:   return(0);
673: }

675: /*@C
676:    NEPSetFunction - Sets the function to compute the nonlinear Function T(lambda)
677:    as well as the location to store the matrix.

679:    Logically Collective on nep

681:    Input Parameters:
682: +  nep - the NEP context
683: .  A   - Function matrix
684: .  B   - preconditioner matrix (usually same as the Function)
685: .  fun - Function evaluation routine (if NULL then NEP retains any
686:          previously set value)
687: -  ctx - [optional] user-defined context for private data for the Function
688:          evaluation routine (may be NULL) (if NULL then NEP retains any
689:          previously set value)

691:    Calling Sequence of fun:
692: $   fun(NEP nep,PetscScalar lambda,Mat F,Mat P,void *ctx)

694: +  nep    - the NEP context
695: .  lambda - the scalar argument where T(.) must be evaluated
696: .  T      - matrix that will contain T(lambda)
697: .  P      - (optional) different matrix to build the preconditioner
698: -  ctx    - (optional) user-defined context, as set by NEPSetFunction()

700:    Level: beginner

702: .seealso: NEPGetFunction(), NEPSetJacobian()
703: @*/
704: PetscErrorCode NEPSetFunction(NEP nep,Mat A,Mat B,PetscErrorCode (*fun)(NEP,PetscScalar,Mat,Mat,void*),void *ctx)
705: {


715:   if (nep->state) { NEPReset(nep); }
716:   else if (nep->fui && nep->fui!=NEP_USER_INTERFACE_CALLBACK) { NEPReset_Problem(nep); }

718:   if (fun) nep->computefunction = fun;
719:   if (ctx) nep->functionctx     = ctx;
720:   if (A) {
721:     PetscObjectReference((PetscObject)A);
722:     MatDestroy(&nep->function);
723:     nep->function = A;
724:   }
725:   if (B) {
726:     PetscObjectReference((PetscObject)B);
727:     MatDestroy(&nep->function_pre);
728:     nep->function_pre = B;
729:   }
730:   nep->fui   = NEP_USER_INTERFACE_CALLBACK;
731:   nep->state = NEP_STATE_INITIAL;
732:   return(0);
733: }

735: /*@C
736:    NEPGetFunction - Returns the Function matrix and optionally the user
737:    provided context for evaluating the Function.

739:    Not Collective, but Mat object will be parallel if NEP object is

741:    Input Parameter:
742: .  nep - the nonlinear eigensolver context

744:    Output Parameters:
745: +  A   - location to stash Function matrix (or NULL)
746: .  B   - location to stash preconditioner matrix (or NULL)
747: .  fun - location to put Function function (or NULL)
748: -  ctx - location to stash Function context (or NULL)

750:    Level: advanced

752: .seealso: NEPSetFunction()
753: @*/
754: PetscErrorCode NEPGetFunction(NEP nep,Mat *A,Mat *B,PetscErrorCode (**fun)(NEP,PetscScalar,Mat,Mat,void*),void **ctx)
755: {
758:   NEPCheckCallback(nep,1);
759:   if (A)   *A   = nep->function;
760:   if (B)   *B   = nep->function_pre;
761:   if (fun) *fun = nep->computefunction;
762:   if (ctx) *ctx = nep->functionctx;
763:   return(0);
764: }

766: /*@C
767:    NEPSetJacobian - Sets the function to compute Jacobian T'(lambda) as well
768:    as the location to store the matrix.

770:    Logically Collective on nep

772:    Input Parameters:
773: +  nep - the NEP context
774: .  A   - Jacobian matrix
775: .  jac - Jacobian evaluation routine (if NULL then NEP retains any
776:          previously set value)
777: -  ctx - [optional] user-defined context for private data for the Jacobian
778:          evaluation routine (may be NULL) (if NULL then NEP retains any
779:          previously set value)

781:    Calling Sequence of jac:
782: $   jac(NEP nep,PetscScalar lambda,Mat J,void *ctx)

784: +  nep    - the NEP context
785: .  lambda - the scalar argument where T'(.) must be evaluated
786: .  J      - matrix that will contain T'(lambda)
787: -  ctx    - (optional) user-defined context, as set by NEPSetJacobian()

789:    Level: beginner

791: .seealso: NEPSetFunction(), NEPGetJacobian()
792: @*/
793: PetscErrorCode NEPSetJacobian(NEP nep,Mat A,PetscErrorCode (*jac)(NEP,PetscScalar,Mat,void*),void *ctx)
794: {


802:   if (nep->state) { NEPReset(nep); }
803:   else if (nep->fui && nep->fui!=NEP_USER_INTERFACE_CALLBACK) { NEPReset_Problem(nep); }

805:   if (jac) nep->computejacobian = jac;
806:   if (ctx) nep->jacobianctx     = ctx;
807:   if (A) {
808:     PetscObjectReference((PetscObject)A);
809:     MatDestroy(&nep->jacobian);
810:     nep->jacobian = A;
811:   }
812:   nep->fui   = NEP_USER_INTERFACE_CALLBACK;
813:   nep->state = NEP_STATE_INITIAL;
814:   return(0);
815: }

817: /*@C
818:    NEPGetJacobian - Returns the Jacobian matrix and optionally the user
819:    provided routine and context for evaluating the Jacobian.

821:    Not Collective, but Mat object will be parallel if NEP object is

823:    Input Parameter:
824: .  nep - the nonlinear eigensolver context

826:    Output Parameters:
827: +  A   - location to stash Jacobian matrix (or NULL)
828: .  jac - location to put Jacobian function (or NULL)
829: -  ctx - location to stash Jacobian context (or NULL)

831:    Level: advanced

833: .seealso: NEPSetJacobian()
834: @*/
835: PetscErrorCode NEPGetJacobian(NEP nep,Mat *A,PetscErrorCode (**jac)(NEP,PetscScalar,Mat,void*),void **ctx)
836: {
839:   NEPCheckCallback(nep,1);
840:   if (A)   *A   = nep->jacobian;
841:   if (jac) *jac = nep->computejacobian;
842:   if (ctx) *ctx = nep->jacobianctx;
843:   return(0);
844: }

846: /*@
847:    NEPSetSplitOperator - Sets the operator of the nonlinear eigenvalue problem
848:    in split form.

850:    Collective on nep

852:    Input Parameters:
853: +  nep - the nonlinear eigensolver context
854: .  n   - number of terms in the split form
855: .  A   - array of matrices
856: .  f   - array of functions
857: -  str - structure flag for matrices

859:    Notes:
860:    The nonlinear operator is written as T(lambda) = sum_i A_i*f_i(lambda),
861:    for i=1,...,n. The derivative T'(lambda) can be obtained using the
862:    derivatives of f_i.

864:    The structure flag provides information about A_i's nonzero pattern
865:    (see MatStructure enum). If all matrices have the same pattern, then
866:    use SAME_NONZERO_PATTERN. If the patterns are different but contained
867:    in the pattern of the first one, then use SUBSET_NONZERO_PATTERN. If
868:    patterns are known to be different, use DIFFERENT_NONZERO_PATTERN.
869:    If set to UNKNOWN_NONZERO_PATTERN, the patterns will be compared to
870:    determine if they are equal.

872:    This function must be called before NEPSetUp(). If it is called again
873:    after NEPSetUp() then the NEP object is reset.

875:    Level: beginner

877: .seealso: NEPGetSplitOperatorTerm(), NEPGetSplitOperatorInfo()
878:  @*/
879: PetscErrorCode NEPSetSplitOperator(NEP nep,PetscInt nt,Mat A[],FN f[],MatStructure str)
880: {
881:   PetscInt       i,n=0,m,m0=0,mloc,nloc,mloc0=0;

887:   if (nt <= 0) SETERRQ1(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more terms, you have %D",nt);

891:   for (i=0;i<nt;i++) {
896:     MatGetSize(A[i],&m,&n);
897:     MatGetLocalSize(A[i],&mloc,&nloc);
898:     if (m!=n) SETERRQ3(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"A[%D] is a non-square matrix (%D rows, %D cols)",i,m,n);
899:     if (mloc!=nloc) SETERRQ3(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"A[%D] does not have equal row and column local sizes (%D, %D)",i,mloc,nloc);
900:     if (!i) { m0 = m; mloc0 = mloc; }
901:     if (m!=m0) SETERRQ3(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_INCOMP,"Dimensions of A[%D] do not match with previous matrices (%D, %D)",i,m,m0);
902:     if (mloc!=mloc0) SETERRQ3(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_INCOMP,"Local dimensions of A[%D] do not match with previous matrices (%D, %D)",i,mloc,mloc0);
903:     PetscObjectReference((PetscObject)A[i]);
904:     PetscObjectReference((PetscObject)f[i]);
905:   }

907:   if (nep->state && (n!=nep->n || nloc!=nep->nloc)) { NEPReset(nep); }
908:   else { NEPReset_Problem(nep); }

910:   /* allocate space and copy matrices and functions */
911:   PetscMalloc1(nt,&nep->A);
912:   PetscLogObjectMemory((PetscObject)nep,nt*sizeof(Mat));
913:   for (i=0;i<nt;i++) nep->A[i] = A[i];
914:   PetscMalloc1(nt,&nep->f);
915:   PetscLogObjectMemory((PetscObject)nep,nt*sizeof(FN));
916:   for (i=0;i<nt;i++) nep->f[i] = f[i];
917:   PetscCalloc1(nt,&nep->nrma);
918:   PetscLogObjectMemory((PetscObject)nep,nt*sizeof(PetscReal));
919:   nep->nt    = nt;
920:   nep->mstr  = str;
921:   nep->fui   = NEP_USER_INTERFACE_SPLIT;
922:   nep->state = NEP_STATE_INITIAL;
923:   return(0);
924: }

926: /*@
927:    NEPGetSplitOperatorTerm - Gets the matrices and functions associated with
928:    the nonlinear operator in split form.

930:    Not collective, though parallel Mats and FNs are returned if the NEP is parallel

932:    Input Parameter:
933: +  nep - the nonlinear eigensolver context
934: -  k   - the index of the requested term (starting in 0)

936:    Output Parameters:
937: +  A - the matrix of the requested term
938: -  f - the function of the requested term

940:    Level: intermediate

942: .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorInfo()
943: @*/
944: PetscErrorCode NEPGetSplitOperatorTerm(NEP nep,PetscInt k,Mat *A,FN *f)
945: {
948:   NEPCheckSplit(nep,1);
949:   if (k<0 || k>=nep->nt) SETERRQ1(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",nep->nt-1);
950:   if (A) *A = nep->A[k];
951:   if (f) *f = nep->f[k];
952:   return(0);
953: }

955: /*@
956:    NEPGetSplitOperatorInfo - Returns the number of terms of the split form of
957:    the nonlinear operator, as well as the structure flag for matrices.

959:    Not collective

961:    Input Parameter:
962: .  nep - the nonlinear eigensolver context

964:    Output Parameters:
965: +  n   - the number of terms passed in NEPSetSplitOperator()
966: -  str - the matrix structure flag passed in NEPSetSplitOperator()

968:    Level: intermediate

970: .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorTerm()
971: @*/
972: PetscErrorCode NEPGetSplitOperatorInfo(NEP nep,PetscInt *n,MatStructure *str)
973: {
976:   NEPCheckSplit(nep,1);
977:   if (n)   *n = nep->nt;
978:   if (str) *str = nep->mstr;
979:   return(0);
980: }