Actual source code: nepbasic.c

slepc-3.14.1 2020-12-08
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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: PetscFunctionList NEPList = 0;
 17: PetscBool         NEPRegisterAllCalled = PETSC_FALSE;
 18: PetscClassId      NEP_CLASSID = 0;
 19: PetscLogEvent     NEP_SetUp = 0,NEP_Solve = 0,NEP_Refine = 0,NEP_FunctionEval = 0,NEP_JacobianEval = 0,NEP_Resolvent = 0;

 21: /*@
 22:    NEPCreate - Creates the default NEP context.

 24:    Collective

 26:    Input Parameter:
 27: .  comm - MPI communicator

 29:    Output Parameter:
 30: .  nep - location to put the NEP context

 32:    Level: beginner

 34: .seealso: NEPSetUp(), NEPSolve(), NEPDestroy(), NEP
 35: @*/
 36: PetscErrorCode NEPCreate(MPI_Comm comm,NEP *outnep)
 37: {
 39:   NEP            nep;

 43:   *outnep = 0;
 44:   NEPInitializePackage();
 45:   SlepcHeaderCreate(nep,NEP_CLASSID,"NEP","Nonlinear Eigenvalue Problem","NEP",comm,NEPDestroy,NEPView);

 47:   nep->max_it          = PETSC_DEFAULT;
 48:   nep->nev             = 1;
 49:   nep->ncv             = PETSC_DEFAULT;
 50:   nep->mpd             = PETSC_DEFAULT;
 51:   nep->nini            = 0;
 52:   nep->target          = 0.0;
 53:   nep->tol             = PETSC_DEFAULT;
 54:   nep->conv            = NEP_CONV_NORM;
 55:   nep->stop            = NEP_STOP_BASIC;
 56:   nep->which           = (NEPWhich)0;
 57:   nep->problem_type    = (NEPProblemType)0;
 58:   nep->refine          = NEP_REFINE_NONE;
 59:   nep->npart           = 1;
 60:   nep->rtol            = PETSC_DEFAULT;
 61:   nep->rits            = PETSC_DEFAULT;
 62:   nep->scheme          = (NEPRefineScheme)0;
 63:   nep->trackall        = PETSC_FALSE;
 64:   nep->twosided        = PETSC_FALSE;

 66:   nep->computefunction = NULL;
 67:   nep->computejacobian = NULL;
 68:   nep->functionctx     = NULL;
 69:   nep->jacobianctx     = NULL;
 70:   nep->converged       = NEPConvergedRelative;
 71:   nep->convergeduser   = NULL;
 72:   nep->convergeddestroy= NULL;
 73:   nep->stopping        = NEPStoppingBasic;
 74:   nep->stoppinguser    = NULL;
 75:   nep->stoppingdestroy = NULL;
 76:   nep->convergedctx    = NULL;
 77:   nep->stoppingctx     = NULL;
 78:   nep->numbermonitors  = 0;

 80:   nep->ds              = NULL;
 81:   nep->V               = NULL;
 82:   nep->W               = NULL;
 83:   nep->rg              = NULL;
 84:   nep->function        = NULL;
 85:   nep->function_pre    = NULL;
 86:   nep->jacobian        = NULL;
 87:   nep->A               = NULL;
 88:   nep->f               = NULL;
 89:   nep->nt              = 0;
 90:   nep->mstr            = DIFFERENT_NONZERO_PATTERN;
 91:   nep->IS              = NULL;
 92:   nep->eigr            = NULL;
 93:   nep->eigi            = NULL;
 94:   nep->errest          = NULL;
 95:   nep->perm            = NULL;
 96:   nep->nwork           = 0;
 97:   nep->work            = NULL;
 98:   nep->data            = NULL;

100:   nep->state           = NEP_STATE_INITIAL;
101:   nep->nconv           = 0;
102:   nep->its             = 0;
103:   nep->n               = 0;
104:   nep->nloc            = 0;
105:   nep->nrma            = NULL;
106:   nep->fui             = (NEPUserInterface)0;
107:   nep->useds           = PETSC_FALSE;
108:   nep->resolvent       = NULL;
109:   nep->reason          = NEP_CONVERGED_ITERATING;

111:   PetscNewLog(nep,&nep->sc);
112:   *outnep = nep;
113:   return(0);
114: }

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

119:    Logically Collective on nep

121:    Input Parameters:
122: +  nep      - the nonlinear eigensolver context
123: -  type     - a known method

125:    Options Database Key:
126: .  -nep_type <method> - Sets the method; use -help for a list
127:     of available methods

129:    Notes:
130:    See "slepc/include/slepcnep.h" for available methods.

132:    Normally, it is best to use the NEPSetFromOptions() command and
133:    then set the NEP type from the options database rather than by using
134:    this routine.  Using the options database provides the user with
135:    maximum flexibility in evaluating the different available methods.
136:    The NEPSetType() routine is provided for those situations where it
137:    is necessary to set the iterative solver independently of the command
138:    line or options database.

140:    Level: intermediate

142: .seealso: NEPType
143: @*/
144: PetscErrorCode NEPSetType(NEP nep,NEPType type)
145: {
146:   PetscErrorCode ierr,(*r)(NEP);
147:   PetscBool      match;


153:   PetscObjectTypeCompare((PetscObject)nep,type,&match);
154:   if (match) return(0);

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

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

162:   nep->state = NEP_STATE_INITIAL;
163:   PetscObjectChangeTypeName((PetscObject)nep,type);
164:   (*r)(nep);
165:   return(0);
166: }

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

171:    Not Collective

173:    Input Parameter:
174: .  nep - the eigensolver context

176:    Output Parameter:
177: .  name - name of NEP method

179:    Level: intermediate

181: .seealso: NEPSetType()
182: @*/
183: PetscErrorCode NEPGetType(NEP nep,NEPType *type)
184: {
188:   *type = ((PetscObject)nep)->type_name;
189:   return(0);
190: }

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

195:    Not Collective

197:    Input Parameters:
198: +  name - name of a new user-defined solver
199: -  function - routine to create the solver context

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

204:    Sample usage:
205: .vb
206:     NEPRegister("my_solver",MySolverCreate);
207: .ve

209:    Then, your solver can be chosen with the procedural interface via
210: $     NEPSetType(nep,"my_solver")
211:    or at runtime via the option
212: $     -nep_type my_solver

214:    Level: advanced

216: .seealso: NEPRegisterAll()
217: @*/
218: PetscErrorCode NEPRegister(const char *name,PetscErrorCode (*function)(NEP))
219: {

223:   NEPInitializePackage();
224:   PetscFunctionListAdd(&NEPList,name,function);
225:   return(0);
226: }

228: /*
229:    NEPReset_Problem - Destroys the problem matrices.
230: @*/
231: PetscErrorCode NEPReset_Problem(NEP nep)
232: {
234:   PetscInt       i;

238:   MatDestroy(&nep->function);
239:   MatDestroy(&nep->function_pre);
240:   MatDestroy(&nep->jacobian);
241:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
242:     MatDestroyMatrices(nep->nt,&nep->A);
243:     for (i=0;i<nep->nt;i++) {
244:       FNDestroy(&nep->f[i]);
245:     }
246:     PetscFree(nep->f);
247:     PetscFree(nep->nrma);
248:     nep->nt = 0;
249:   }
250:   return(0);
251: }
252: /*@
253:    NEPReset - Resets the NEP context to the initial state (prior to setup)
254:    and destroys any allocated Vecs and Mats.

256:    Collective on nep

258:    Input Parameter:
259: .  nep - eigensolver context obtained from NEPCreate()

261:    Level: advanced

263: .seealso: NEPDestroy()
264: @*/
265: PetscErrorCode NEPReset(NEP nep)
266: {

271:   if (!nep) return(0);
272:   if (nep->ops->reset) { (nep->ops->reset)(nep); }
273:   if (nep->refineksp) { KSPReset(nep->refineksp); }
274:   NEPReset_Problem(nep);
275:   BVDestroy(&nep->V);
276:   BVDestroy(&nep->W);
277:   VecDestroyVecs(nep->nwork,&nep->work);
278:   MatDestroy(&nep->resolvent);
279:   nep->nwork = 0;
280:   nep->state = NEP_STATE_INITIAL;
281:   return(0);
282: }

284: /*@
285:    NEPDestroy - Destroys the NEP context.

287:    Collective on nep

289:    Input Parameter:
290: .  nep - eigensolver context obtained from NEPCreate()

292:    Level: beginner

294: .seealso: NEPCreate(), NEPSetUp(), NEPSolve()
295: @*/
296: PetscErrorCode NEPDestroy(NEP *nep)
297: {

301:   if (!*nep) return(0);
303:   if (--((PetscObject)(*nep))->refct > 0) { *nep = 0; return(0); }
304:   NEPReset(*nep);
305:   if ((*nep)->ops->destroy) { (*(*nep)->ops->destroy)(*nep); }
306:   if ((*nep)->eigr) {
307:     PetscFree4((*nep)->eigr,(*nep)->eigi,(*nep)->errest,(*nep)->perm);
308:   }
309:   RGDestroy(&(*nep)->rg);
310:   DSDestroy(&(*nep)->ds);
311:   KSPDestroy(&(*nep)->refineksp);
312:   PetscSubcommDestroy(&(*nep)->refinesubc);
313:   PetscFree((*nep)->sc);
314:   /* just in case the initial vectors have not been used */
315:   SlepcBasisDestroy_Private(&(*nep)->nini,&(*nep)->IS);
316:   if ((*nep)->convergeddestroy) {
317:     (*(*nep)->convergeddestroy)((*nep)->convergedctx);
318:   }
319:   NEPMonitorCancel(*nep);
320:   PetscHeaderDestroy(nep);
321:   return(0);
322: }

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

327:    Collective on nep

329:    Input Parameters:
330: +  nep - eigensolver context obtained from NEPCreate()
331: -  bv  - the basis vectors object

333:    Note:
334:    Use NEPGetBV() to retrieve the basis vectors context (for example,
335:    to free it at the end of the computations).

337:    Level: advanced

339: .seealso: NEPGetBV()
340: @*/
341: PetscErrorCode NEPSetBV(NEP nep,BV bv)
342: {

349:   PetscObjectReference((PetscObject)bv);
350:   BVDestroy(&nep->V);
351:   nep->V = bv;
352:   PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->V);
353:   return(0);
354: }

356: /*@
357:    NEPGetBV - Obtain the basis vectors object associated to the nonlinear
358:    eigensolver object.

360:    Not Collective

362:    Input Parameters:
363: .  nep - eigensolver context obtained from NEPCreate()

365:    Output Parameter:
366: .  bv - basis vectors context

368:    Level: advanced

370: .seealso: NEPSetBV()
371: @*/
372: PetscErrorCode NEPGetBV(NEP nep,BV *bv)
373: {

379:   if (!nep->V) {
380:     BVCreate(PetscObjectComm((PetscObject)nep),&nep->V);
381:     PetscObjectIncrementTabLevel((PetscObject)nep->V,(PetscObject)nep,0);
382:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->V);
383:     PetscObjectSetOptions((PetscObject)nep->V,((PetscObject)nep)->options);
384:   }
385:   *bv = nep->V;
386:   return(0);
387: }

389: /*@
390:    NEPSetRG - Associates a region object to the nonlinear eigensolver.

392:    Collective on nep

394:    Input Parameters:
395: +  nep - eigensolver context obtained from NEPCreate()
396: -  rg  - the region object

398:    Note:
399:    Use NEPGetRG() to retrieve the region context (for example,
400:    to free it at the end of the computations).

402:    Level: advanced

404: .seealso: NEPGetRG()
405: @*/
406: PetscErrorCode NEPSetRG(NEP nep,RG rg)
407: {

414:   PetscObjectReference((PetscObject)rg);
415:   RGDestroy(&nep->rg);
416:   nep->rg = rg;
417:   PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->rg);
418:   return(0);
419: }

421: /*@
422:    NEPGetRG - Obtain the region object associated to the
423:    nonlinear eigensolver object.

425:    Not Collective

427:    Input Parameters:
428: .  nep - eigensolver context obtained from NEPCreate()

430:    Output Parameter:
431: .  rg - region context

433:    Level: advanced

435: .seealso: NEPSetRG()
436: @*/
437: PetscErrorCode NEPGetRG(NEP nep,RG *rg)
438: {

444:   if (!nep->rg) {
445:     RGCreate(PetscObjectComm((PetscObject)nep),&nep->rg);
446:     PetscObjectIncrementTabLevel((PetscObject)nep->rg,(PetscObject)nep,0);
447:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->rg);
448:     PetscObjectSetOptions((PetscObject)nep->rg,((PetscObject)nep)->options);
449:   }
450:   *rg = nep->rg;
451:   return(0);
452: }

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

457:    Collective on nep

459:    Input Parameters:
460: +  nep - eigensolver context obtained from NEPCreate()
461: -  ds  - the direct solver object

463:    Note:
464:    Use NEPGetDS() to retrieve the direct solver context (for example,
465:    to free it at the end of the computations).

467:    Level: advanced

469: .seealso: NEPGetDS()
470: @*/
471: PetscErrorCode NEPSetDS(NEP nep,DS ds)
472: {

479:   PetscObjectReference((PetscObject)ds);
480:   DSDestroy(&nep->ds);
481:   nep->ds = ds;
482:   PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->ds);
483:   return(0);
484: }

486: /*@
487:    NEPGetDS - Obtain the direct solver object associated to the
488:    nonlinear eigensolver object.

490:    Not Collective

492:    Input Parameters:
493: .  nep - eigensolver context obtained from NEPCreate()

495:    Output Parameter:
496: .  ds - direct solver context

498:    Level: advanced

500: .seealso: NEPSetDS()
501: @*/
502: PetscErrorCode NEPGetDS(NEP nep,DS *ds)
503: {

509:   if (!nep->ds) {
510:     DSCreate(PetscObjectComm((PetscObject)nep),&nep->ds);
511:     PetscObjectIncrementTabLevel((PetscObject)nep->ds,(PetscObject)nep,0);
512:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->ds);
513:     PetscObjectSetOptions((PetscObject)nep->ds,((PetscObject)nep)->options);
514:   }
515:   *ds = nep->ds;
516:   return(0);
517: }

519: /*@
520:    NEPRefineGetKSP - Obtain the ksp object used by the eigensolver
521:    object in the refinement phase.

523:    Not Collective

525:    Input Parameters:
526: .  nep - eigensolver context obtained from NEPCreate()

528:    Output Parameter:
529: .  ksp - ksp context

531:    Level: advanced

533: .seealso: NEPSetRefine()
534: @*/
535: PetscErrorCode NEPRefineGetKSP(NEP nep,KSP *ksp)
536: {

542:   if (!nep->refineksp) {
543:     if (nep->npart>1) {
544:       /* Split in subcomunicators */
545:       PetscSubcommCreate(PetscObjectComm((PetscObject)nep),&nep->refinesubc);
546:       PetscSubcommSetNumber(nep->refinesubc,nep->npart);
547:       PetscSubcommSetType(nep->refinesubc,PETSC_SUBCOMM_CONTIGUOUS);
548:       PetscLogObjectMemory((PetscObject)nep,sizeof(PetscSubcomm));
549:     }
550:     KSPCreate((nep->npart==1)?PetscObjectComm((PetscObject)nep):PetscSubcommChild(nep->refinesubc),&nep->refineksp);
551:     PetscObjectIncrementTabLevel((PetscObject)nep->refineksp,(PetscObject)nep,0);
552:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->refineksp);
553:     PetscObjectSetOptions((PetscObject)nep->refineksp,((PetscObject)nep)->options);
554:     KSPSetOptionsPrefix(*ksp,((PetscObject)nep)->prefix);
555:     KSPAppendOptionsPrefix(*ksp,"nep_refine_");
556:     KSPSetTolerances(nep->refineksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
557:   }
558:   *ksp = nep->refineksp;
559:   return(0);
560: }

562: /*@
563:    NEPSetTarget - Sets the value of the target.

565:    Logically Collective on nep

567:    Input Parameters:
568: +  nep    - eigensolver context
569: -  target - the value of the target

571:    Options Database Key:
572: .  -nep_target <scalar> - the value of the target

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

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

582:    Level: intermediate

584: .seealso: NEPGetTarget(), NEPSetWhichEigenpairs()
585: @*/
586: PetscErrorCode NEPSetTarget(NEP nep,PetscScalar target)
587: {
591:   nep->target = target;
592:   return(0);
593: }

595: /*@
596:    NEPGetTarget - Gets the value of the target.

598:    Not Collective

600:    Input Parameter:
601: .  nep - eigensolver context

603:    Output Parameter:
604: .  target - the value of the target

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

609:    Level: intermediate

611: .seealso: NEPSetTarget()
612: @*/
613: PetscErrorCode NEPGetTarget(NEP nep,PetscScalar* target)
614: {
618:   *target = nep->target;
619:   return(0);
620: }

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

626:    Logically Collective on nep

628:    Input Parameters:
629: +  nep - the NEP context
630: .  A   - Function matrix
631: .  B   - preconditioner matrix (usually same as the Function)
632: .  fun - Function evaluation routine (if NULL then NEP retains any
633:          previously set value)
634: -  ctx - [optional] user-defined context for private data for the Function
635:          evaluation routine (may be NULL) (if NULL then NEP retains any
636:          previously set value)

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

641: +  nep    - the NEP context
642: .  lambda - the scalar argument where T(.) must be evaluated
643: .  T      - matrix that will contain T(lambda)
644: .  P      - (optional) different matrix to build the preconditioner
645: -  ctx    - (optional) user-defined context, as set by NEPSetFunction()

647:    Level: beginner

649: .seealso: NEPGetFunction(), NEPSetJacobian()
650: @*/
651: PetscErrorCode NEPSetFunction(NEP nep,Mat A,Mat B,PetscErrorCode (*fun)(NEP,PetscScalar,Mat,Mat,void*),void *ctx)
652: {


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

665:   if (fun) nep->computefunction = fun;
666:   if (ctx) nep->functionctx     = ctx;
667:   if (A) {
668:     PetscObjectReference((PetscObject)A);
669:     MatDestroy(&nep->function);
670:     nep->function = A;
671:   }
672:   if (B) {
673:     PetscObjectReference((PetscObject)B);
674:     MatDestroy(&nep->function_pre);
675:     nep->function_pre = B;
676:   }
677:   nep->fui   = NEP_USER_INTERFACE_CALLBACK;
678:   nep->state = NEP_STATE_INITIAL;
679:   return(0);
680: }

682: /*@C
683:    NEPGetFunction - Returns the Function matrix and optionally the user
684:    provided context for evaluating the Function.

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

688:    Input Parameter:
689: .  nep - the nonlinear eigensolver context

691:    Output Parameters:
692: +  A   - location to stash Function matrix (or NULL)
693: .  B   - location to stash preconditioner matrix (or NULL)
694: .  fun - location to put Function function (or NULL)
695: -  ctx - location to stash Function context (or NULL)

697:    Level: advanced

699: .seealso: NEPSetFunction()
700: @*/
701: PetscErrorCode NEPGetFunction(NEP nep,Mat *A,Mat *B,PetscErrorCode (**fun)(NEP,PetscScalar,Mat,Mat,void*),void **ctx)
702: {
705:   NEPCheckCallback(nep,1);
706:   if (A)   *A   = nep->function;
707:   if (B)   *B   = nep->function_pre;
708:   if (fun) *fun = nep->computefunction;
709:   if (ctx) *ctx = nep->functionctx;
710:   return(0);
711: }

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

717:    Logically Collective on nep

719:    Input Parameters:
720: +  nep - the NEP context
721: .  A   - Jacobian matrix
722: .  jac - Jacobian evaluation routine (if NULL then NEP retains any
723:          previously set value)
724: -  ctx - [optional] user-defined context for private data for the Jacobian
725:          evaluation routine (may be NULL) (if NULL then NEP retains any
726:          previously set value)

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

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

736:    Level: beginner

738: .seealso: NEPSetFunction(), NEPGetJacobian()
739: @*/
740: PetscErrorCode NEPSetJacobian(NEP nep,Mat A,PetscErrorCode (*jac)(NEP,PetscScalar,Mat,void*),void *ctx)
741: {


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

752:   if (jac) nep->computejacobian = jac;
753:   if (ctx) nep->jacobianctx     = ctx;
754:   if (A) {
755:     PetscObjectReference((PetscObject)A);
756:     MatDestroy(&nep->jacobian);
757:     nep->jacobian = A;
758:   }
759:   nep->fui   = NEP_USER_INTERFACE_CALLBACK;
760:   nep->state = NEP_STATE_INITIAL;
761:   return(0);
762: }

764: /*@C
765:    NEPGetJacobian - Returns the Jacobian matrix and optionally the user
766:    provided routine and context for evaluating the Jacobian.

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

770:    Input Parameter:
771: .  nep - the nonlinear eigensolver context

773:    Output Parameters:
774: +  A   - location to stash Jacobian matrix (or NULL)
775: .  jac - location to put Jacobian function (or NULL)
776: -  ctx - location to stash Jacobian context (or NULL)

778:    Level: advanced

780: .seealso: NEPSetJacobian()
781: @*/
782: PetscErrorCode NEPGetJacobian(NEP nep,Mat *A,PetscErrorCode (**jac)(NEP,PetscScalar,Mat,void*),void **ctx)
783: {
786:   NEPCheckCallback(nep,1);
787:   if (A)   *A   = nep->jacobian;
788:   if (jac) *jac = nep->computejacobian;
789:   if (ctx) *ctx = nep->jacobianctx;
790:   return(0);
791: }

793: /*@
794:    NEPSetSplitOperator - Sets the operator of the nonlinear eigenvalue problem
795:    in split form.

797:    Collective on nep

799:    Input Parameters:
800: +  nep - the nonlinear eigensolver context
801: .  n   - number of terms in the split form
802: .  A   - array of matrices
803: .  f   - array of functions
804: -  str - structure flag for matrices

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

811:    The structure flag provides information about A_i's nonzero pattern
812:    (see MatStructure enum). If all matrices have the same pattern, then
813:    use SAME_NONZERO_PATTERN. If the patterns are different but contained
814:    in the pattern of the first one, then use SUBSET_NONZERO_PATTERN.
815:    Otherwise use DIFFERENT_NONZERO_PATTERN.

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

820:    Level: beginner

822: .seealso: NEPGetSplitOperatorTerm(), NEPGetSplitOperatorInfo()
823:  @*/
824: PetscErrorCode NEPSetSplitOperator(NEP nep,PetscInt n,Mat A[],FN f[],MatStructure str)
825: {
826:   PetscInt       i;

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

838:   for (i=0;i<n;i++) {
841:     PetscObjectReference((PetscObject)A[i]);
842:     PetscObjectReference((PetscObject)f[i]);
843:   }

845:   if (nep->state) { NEPReset(nep); }
846:   else { NEPReset_Problem(nep); }

848:   /* allocate space and copy matrices and functions */
849:   PetscMalloc1(n,&nep->A);
850:   PetscLogObjectMemory((PetscObject)nep,n*sizeof(Mat));
851:   for (i=0;i<n;i++) nep->A[i] = A[i];
852:   PetscMalloc1(n,&nep->f);
853:   PetscLogObjectMemory((PetscObject)nep,n*sizeof(FN));
854:   for (i=0;i<n;i++) nep->f[i] = f[i];
855:   PetscCalloc1(n,&nep->nrma);
856:   PetscLogObjectMemory((PetscObject)nep,n*sizeof(PetscReal));
857:   nep->nt    = n;
858:   nep->mstr  = str;
859:   nep->fui   = NEP_USER_INTERFACE_SPLIT;
860:   nep->state = NEP_STATE_INITIAL;
861:   return(0);
862: }

864: /*@
865:    NEPGetSplitOperatorTerm - Gets the matrices and functions associated with
866:    the nonlinear operator in split form.

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

870:    Input Parameter:
871: +  nep - the nonlinear eigensolver context
872: -  k   - the index of the requested term (starting in 0)

874:    Output Parameters:
875: +  A - the matrix of the requested term
876: -  f - the function of the requested term

878:    Level: intermediate

880: .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorInfo()
881: @*/
882: PetscErrorCode NEPGetSplitOperatorTerm(NEP nep,PetscInt k,Mat *A,FN *f)
883: {
886:   NEPCheckSplit(nep,1);
887:   if (k<0 || k>=nep->nt) SETERRQ1(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",nep->nt-1);
888:   if (A) *A = nep->A[k];
889:   if (f) *f = nep->f[k];
890:   return(0);
891: }

893: /*@
894:    NEPGetSplitOperatorInfo - Returns the number of terms of the split form of
895:    the nonlinear operator, as well as the structure flag for matrices.

897:    Not collective

899:    Input Parameter:
900: .  nep - the nonlinear eigensolver context

902:    Output Parameters:
903: +  n   - the number of terms passed in NEPSetSplitOperator()
904: -  str - the matrix structure flag passed in NEPSetSplitOperator()

906:    Level: intermediate

908: .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorTerm()
909: @*/
910: PetscErrorCode NEPGetSplitOperatorInfo(NEP nep,PetscInt *n,MatStructure *str)
911: {
914:   NEPCheckSplit(nep,1);
915:   if (n)   *n = nep->nt;
916:   if (str) *str = nep->mstr;
917:   return(0);
918: }