1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, 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> /*I "slepcnep.h" I*/
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_DerivativesEval = 0;
21: /*@
22: NEPCreate - Creates the default NEP context.
24: Collective on MPI_Comm
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 = 0;
48: nep->nev = 1;
49: nep->ncv = 0;
50: nep->mpd = 0;
51: nep->nini = 0;
52: nep->target = 0.0;
53: nep->tol = PETSC_DEFAULT;
54: nep->conv = NEP_CONV_REL;
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;
65: nep->computefunction = NULL;
66: nep->computejacobian = NULL;
67: nep->functionctx = NULL;
68: nep->jacobianctx = NULL;
69: nep->computederivatives = NULL;
70: nep->derivativesctx = NULL;
71: nep->converged = NEPConvergedRelative;
72: nep->convergeduser = NULL;
73: nep->convergeddestroy= NULL;
74: nep->stopping = NEPStoppingBasic;
75: nep->stoppinguser = NULL;
76: nep->stoppingdestroy = NULL;
77: nep->convergedctx = NULL;
78: nep->stoppingctx = NULL;
79: nep->numbermonitors = 0;
81: nep->ds = NULL;
82: nep->V = NULL;
83: nep->rg = NULL;
84: nep->function = NULL;
85: nep->function_pre = NULL;
86: nep->jacobian = NULL;
87: nep->derivatives = NULL;
88: nep->A = NULL;
89: nep->f = NULL;
90: nep->nt = 0;
91: nep->mstr = DIFFERENT_NONZERO_PATTERN;
92: nep->IS = NULL;
93: nep->eigr = NULL;
94: nep->eigi = NULL;
95: nep->errest = NULL;
96: nep->perm = NULL;
97: nep->nwork = 0;
98: nep->work = NULL;
99: nep->data = NULL;
101: nep->state = NEP_STATE_INITIAL;
102: nep->nconv = 0;
103: nep->its = 0;
104: nep->n = 0;
105: nep->nloc = 0;
106: nep->nrma = NULL;
107: nep->fui = (NEPUserInterface)0;
108: nep->reason = NEP_CONVERGED_ITERATING;
110: PetscNewLog(nep,&nep->sc);
111: *outnep = nep;
112: return(0);
113: }
115: /*@C
116: NEPSetType - Selects the particular solver to be used in the NEP object.
118: Logically Collective on NEP120: Input Parameters:
121: + nep - the nonlinear eigensolver context
122: - type - a known method
124: Options Database Key:
125: . -nep_type <method> - Sets the method; use -help for a list
126: of available methods
128: Notes:
129: See "slepc/include/slepcnep.h" for available methods.
131: Normally, it is best to use the NEPSetFromOptions() command and
132: then set the NEP type from the options database rather than by using
133: this routine. Using the options database provides the user with
134: maximum flexibility in evaluating the different available methods.
135: The NEPSetType() routine is provided for those situations where it
136: is necessary to set the iterative solver independently of the command
137: line or options database.
139: Level: intermediate
141: .seealso: NEPType142: @*/
143: PetscErrorCode NEPSetType(NEP nep,NEPType type)144: {
145: PetscErrorCode ierr,(*r)(NEP);
146: PetscBool match;
152: PetscObjectTypeCompare((PetscObject)nep,type,&match);
153: if (match) return(0);
155: PetscFunctionListFind(NEPList,type,&r);
156: if (!r) SETERRQ1(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown NEP type given: %s",type);
158: if (nep->ops->destroy) { (*nep->ops->destroy)(nep); }
159: PetscMemzero(nep->ops,sizeof(struct _NEPOps));
161: nep->state = NEP_STATE_INITIAL;
162: PetscObjectChangeTypeName((PetscObject)nep,type);
163: (*r)(nep);
164: return(0);
165: }
167: /*@C
168: NEPGetType - Gets the NEP type as a string from the NEP object.
170: Not Collective
172: Input Parameter:
173: . nep - the eigensolver context
175: Output Parameter:
176: . name - name of NEP method
178: Level: intermediate
180: .seealso: NEPSetType()
181: @*/
182: PetscErrorCode NEPGetType(NEP nep,NEPType *type)183: {
187: *type = ((PetscObject)nep)->type_name;
188: return(0);
189: }
191: /*@C
192: NEPRegister - Adds a method to the nonlinear eigenproblem solver package.
194: Not Collective
196: Input Parameters:
197: + name - name of a new user-defined solver
198: - function - routine to create the solver context
200: Notes:
201: NEPRegister() may be called multiple times to add several user-defined solvers.
203: Sample usage:
204: .vb
205: NEPRegister("my_solver",MySolverCreate);
206: .ve
208: Then, your solver can be chosen with the procedural interface via
209: $ NEPSetType(nep,"my_solver")
210: or at runtime via the option
211: $ -nep_type my_solver
213: Level: advanced
215: .seealso: NEPRegisterAll()
216: @*/
217: PetscErrorCode NEPRegister(const char *name,PetscErrorCode (*function)(NEP))218: {
222: PetscFunctionListAdd(&NEPList,name,function);
223: return(0);
224: }
226: /*
227: NEPReset_Problem - Destroys the problem matrices.
228: @*/
229: PetscErrorCode NEPReset_Problem(NEP nep)230: {
232: PetscInt i;
236: MatDestroy(&nep->function);
237: MatDestroy(&nep->function_pre);
238: MatDestroy(&nep->jacobian);
239: MatDestroy(&nep->derivatives);
240: if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
241: MatDestroyMatrices(nep->nt,&nep->A);
242: for (i=0;i<nep->nt;i++) {
243: FNDestroy(&nep->f[i]);
244: }
245: PetscFree(nep->f);
246: PetscFree(nep->nrma);
247: nep->nt = 0;
248: }
249: return(0);
250: }
251: /*@
252: NEPReset - Resets the NEP context to the initial state (prior to setup)
253: and destroys any allocated Vecs and Mats.
255: Collective on NEP257: Input Parameter:
258: . nep - eigensolver context obtained from NEPCreate()
260: Level: advanced
262: .seealso: NEPDestroy()
263: @*/
264: PetscErrorCode NEPReset(NEP nep)265: {
270: if (!nep) return(0);
271: if (nep->ops->reset) { (nep->ops->reset)(nep); }
272: if (nep->refineksp) { KSPReset(nep->refineksp); }
273: NEPReset_Problem(nep);
274: BVDestroy(&nep->V);
275: VecDestroyVecs(nep->nwork,&nep->work);
276: nep->nwork = 0;
277: nep->state = NEP_STATE_INITIAL;
278: return(0);
279: }
281: /*@
282: NEPDestroy - Destroys the NEP context.
284: Collective on NEP286: Input Parameter:
287: . nep - eigensolver context obtained from NEPCreate()
289: Level: beginner
291: .seealso: NEPCreate(), NEPSetUp(), NEPSolve()
292: @*/
293: PetscErrorCode NEPDestroy(NEP *nep)294: {
298: if (!*nep) return(0);
300: if (--((PetscObject)(*nep))->refct > 0) { *nep = 0; return(0); }
301: NEPReset(*nep);
302: if ((*nep)->ops->destroy) { (*(*nep)->ops->destroy)(*nep); }
303: if ((*nep)->eigr) {
304: PetscFree4((*nep)->eigr,(*nep)->eigi,(*nep)->errest,(*nep)->perm);
305: }
306: RGDestroy(&(*nep)->rg);
307: DSDestroy(&(*nep)->ds);
308: KSPDestroy(&(*nep)->refineksp);
309: PetscSubcommDestroy(&(*nep)->refinesubc);
310: PetscFree((*nep)->sc);
311: /* just in case the initial vectors have not been used */
312: SlepcBasisDestroy_Private(&(*nep)->nini,&(*nep)->IS);
313: if ((*nep)->convergeddestroy) {
314: (*(*nep)->convergeddestroy)((*nep)->convergedctx);
315: }
316: NEPMonitorCancel(*nep);
317: PetscHeaderDestroy(nep);
318: return(0);
319: }
321: /*@
322: NEPSetBV - Associates a basis vectors object to the nonlinear eigensolver.
324: Collective on NEP326: Input Parameters:
327: + nep - eigensolver context obtained from NEPCreate()
328: - bv - the basis vectors object
330: Note:
331: Use NEPGetBV() to retrieve the basis vectors context (for example,
332: to free it at the end of the computations).
334: Level: advanced
336: .seealso: NEPGetBV()
337: @*/
338: PetscErrorCode NEPSetBV(NEP nep,BV bv)339: {
346: PetscObjectReference((PetscObject)bv);
347: BVDestroy(&nep->V);
348: nep->V = bv;
349: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->V);
350: return(0);
351: }
353: /*@
354: NEPGetBV - Obtain the basis vectors object associated to the nonlinear
355: eigensolver object.
357: Not Collective
359: Input Parameters:
360: . nep - eigensolver context obtained from NEPCreate()
362: Output Parameter:
363: . bv - basis vectors context
365: Level: advanced
367: .seealso: NEPSetBV()
368: @*/
369: PetscErrorCode NEPGetBV(NEP nep,BV *bv)370: {
376: if (!nep->V) {
377: BVCreate(PetscObjectComm((PetscObject)nep),&nep->V);
378: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->V);
379: }
380: *bv = nep->V;
381: return(0);
382: }
384: /*@
385: NEPSetRG - Associates a region object to the nonlinear eigensolver.
387: Collective on NEP389: Input Parameters:
390: + nep - eigensolver context obtained from NEPCreate()
391: - rg - the region object
393: Note:
394: Use NEPGetRG() to retrieve the region context (for example,
395: to free it at the end of the computations).
397: Level: advanced
399: .seealso: NEPGetRG()
400: @*/
401: PetscErrorCode NEPSetRG(NEP nep,RG rg)402: {
409: PetscObjectReference((PetscObject)rg);
410: RGDestroy(&nep->rg);
411: nep->rg = rg;
412: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->rg);
413: return(0);
414: }
416: /*@
417: NEPGetRG - Obtain the region object associated to the
418: nonlinear eigensolver object.
420: Not Collective
422: Input Parameters:
423: . nep - eigensolver context obtained from NEPCreate()
425: Output Parameter:
426: . rg - region context
428: Level: advanced
430: .seealso: NEPSetRG()
431: @*/
432: PetscErrorCode NEPGetRG(NEP nep,RG *rg)433: {
439: if (!nep->rg) {
440: RGCreate(PetscObjectComm((PetscObject)nep),&nep->rg);
441: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->rg);
442: }
443: *rg = nep->rg;
444: return(0);
445: }
447: /*@
448: NEPSetDS - Associates a direct solver object to the nonlinear eigensolver.
450: Collective on NEP452: Input Parameters:
453: + nep - eigensolver context obtained from NEPCreate()
454: - ds - the direct solver object
456: Note:
457: Use NEPGetDS() to retrieve the direct solver context (for example,
458: to free it at the end of the computations).
460: Level: advanced
462: .seealso: NEPGetDS()
463: @*/
464: PetscErrorCode NEPSetDS(NEP nep,DS ds)465: {
472: PetscObjectReference((PetscObject)ds);
473: DSDestroy(&nep->ds);
474: nep->ds = ds;
475: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->ds);
476: return(0);
477: }
479: /*@
480: NEPGetDS - Obtain the direct solver object associated to the
481: nonlinear eigensolver object.
483: Not Collective
485: Input Parameters:
486: . nep - eigensolver context obtained from NEPCreate()
488: Output Parameter:
489: . ds - direct solver context
491: Level: advanced
493: .seealso: NEPSetDS()
494: @*/
495: PetscErrorCode NEPGetDS(NEP nep,DS *ds)496: {
502: if (!nep->ds) {
503: DSCreate(PetscObjectComm((PetscObject)nep),&nep->ds);
504: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->ds);
505: }
506: *ds = nep->ds;
507: return(0);
508: }
510: /*@
511: NEPRefineGetKSP - Obtain the ksp object used by the eigensolver
512: object in the refinement phase.
514: Not Collective
516: Input Parameters:
517: . nep - eigensolver context obtained from NEPCreate()
519: Output Parameter:
520: . ksp - ksp context
522: Level: advanced
524: .seealso: NEPSetRefine()
525: @*/
526: PetscErrorCode NEPRefineGetKSP(NEP nep,KSP *ksp)527: {
533: if (!nep->refineksp) {
534: if (nep->npart>1) {
535: /* Split in subcomunicators */
536: PetscSubcommCreate(PetscObjectComm((PetscObject)nep),&nep->refinesubc);
537: PetscSubcommSetNumber(nep->refinesubc,nep->npart);
538: PetscSubcommSetType(nep->refinesubc,PETSC_SUBCOMM_CONTIGUOUS);
539: PetscLogObjectMemory((PetscObject)nep,sizeof(PetscSubcomm));
540: }
541: KSPCreate((nep->npart==1)?PetscObjectComm((PetscObject)nep):PetscSubcommChild(nep->refinesubc),&nep->refineksp);
542: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->refineksp);
543: KSPSetOptionsPrefix(*ksp,((PetscObject)nep)->prefix);
544: KSPAppendOptionsPrefix(*ksp,"nep_refine_");
545: KSPSetTolerances(nep->refineksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
546: }
547: *ksp = nep->refineksp;
548: return(0);
549: }
551: /*@
552: NEPSetTarget - Sets the value of the target.
554: Logically Collective on NEP556: Input Parameters:
557: + nep - eigensolver context
558: - target - the value of the target
560: Options Database Key:
561: . -nep_target <scalar> - the value of the target
563: Notes:
564: The target is a scalar value used to determine the portion of the spectrum
565: of interest. It is used in combination with NEPSetWhichEigenpairs().
567: In the case of complex scalars, a complex value can be provided in the
568: command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
569: -nep_target 1.0+2.0i
571: Level: intermediate
573: .seealso: NEPGetTarget(), NEPSetWhichEigenpairs()
574: @*/
575: PetscErrorCode NEPSetTarget(NEP nep,PetscScalar target)576: {
580: nep->target = target;
581: return(0);
582: }
584: /*@
585: NEPGetTarget - Gets the value of the target.
587: Not Collective
589: Input Parameter:
590: . nep - eigensolver context
592: Output Parameter:
593: . target - the value of the target
595: Note:
596: If the target was not set by the user, then zero is returned.
598: Level: intermediate
600: .seealso: NEPSetTarget()
601: @*/
602: PetscErrorCode NEPGetTarget(NEP nep,PetscScalar* target)603: {
607: *target = nep->target;
608: return(0);
609: }
611: /*@C
612: NEPSetFunction - Sets the function to compute the nonlinear Function T(lambda)
613: as well as the location to store the matrix.
615: Logically Collective on NEP and Mat
617: Input Parameters:
618: + nep - the NEP context
619: . A - Function matrix
620: . B - preconditioner matrix (usually same as the Function)
621: . fun - Function evaluation routine (if NULL then NEP retains any
622: previously set value)
623: - ctx - [optional] user-defined context for private data for the Function
624: evaluation routine (may be NULL) (if NULL then NEP retains any
625: previously set value)
627: Calling Sequence of fun:
628: $ fun(NEP nep,PetscScalar lambda,Mat F,Mat P,void *ctx)
630: + nep - the NEP context
631: . lambda - the scalar argument where T(.) must be evaluated
632: . T - matrix that will contain T(lambda)
633: . P - (optional) different matrix to build the preconditioner
634: - ctx - (optional) user-defined context, as set by NEPSetFunction()
636: Level: beginner
638: .seealso: NEPGetFunction(), NEPSetJacobian()
639: @*/
640: PetscErrorCode NEPSetFunction(NEP nep,Mat A,Mat B,PetscErrorCode (*fun)(NEP,PetscScalar,Mat,Mat,void*),void *ctx)641: {
651: if (nep->state) { NEPReset(nep); }
652: else if (nep->fui && nep->fui!=NEP_USER_INTERFACE_CALLBACK) { NEPReset_Problem(nep); }
654: if (fun) nep->computefunction = fun;
655: if (ctx) nep->functionctx = ctx;
656: if (A) {
657: PetscObjectReference((PetscObject)A);
658: MatDestroy(&nep->function);
659: nep->function = A;
660: }
661: if (B) {
662: PetscObjectReference((PetscObject)B);
663: MatDestroy(&nep->function_pre);
664: nep->function_pre = B;
665: }
666: nep->fui = NEP_USER_INTERFACE_CALLBACK;
667: nep->state = NEP_STATE_INITIAL;
668: return(0);
669: }
671: /*@C
672: NEPGetFunction - Returns the Function matrix and optionally the user
673: provided context for evaluating the Function.
675: Not Collective, but Mat object will be parallel if NEP object is
677: Input Parameter:
678: . nep - the nonlinear eigensolver context
680: Output Parameters:
681: + A - location to stash Function matrix (or NULL)
682: . B - location to stash preconditioner matrix (or NULL)
683: . fun - location to put Function function (or NULL)
684: - ctx - location to stash Function context (or NULL)
686: Level: advanced
688: .seealso: NEPSetFunction()
689: @*/
690: PetscErrorCode NEPGetFunction(NEP nep,Mat *A,Mat *B,PetscErrorCode (**fun)(NEP,PetscScalar,Mat,Mat,void*),void **ctx)691: {
694: NEPCheckCallback(nep,1);
695: if (A) *A = nep->function;
696: if (B) *B = nep->function_pre;
697: if (fun) *fun = nep->computefunction;
698: if (ctx) *ctx = nep->functionctx;
699: return(0);
700: }
702: /*@C
703: NEPSetJacobian - Sets the function to compute Jacobian T'(lambda) as well
704: as the location to store the matrix.
706: Logically Collective on NEP and Mat
708: Input Parameters:
709: + nep - the NEP context
710: . A - Jacobian matrix
711: . jac - Jacobian evaluation routine (if NULL then NEP retains any
712: previously set value)
713: - ctx - [optional] user-defined context for private data for the Jacobian
714: evaluation routine (may be NULL) (if NULL then NEP retains any
715: previously set value)
717: Calling Sequence of jac:
718: $ jac(NEP nep,PetscScalar lambda,Mat J,void *ctx)
720: + nep - the NEP context
721: . lambda - the scalar argument where T'(.) must be evaluated
722: . J - matrix that will contain T'(lambda)
723: - ctx - (optional) user-defined context, as set by NEPSetJacobian()
725: Level: beginner
727: .seealso: NEPSetFunction(), NEPGetJacobian()
728: @*/
729: PetscErrorCode NEPSetJacobian(NEP nep,Mat A,PetscErrorCode (*jac)(NEP,PetscScalar,Mat,void*),void *ctx)730: {
738: if (nep->state) { NEPReset(nep); }
739: else if (nep->fui && nep->fui!=NEP_USER_INTERFACE_CALLBACK) { NEPReset_Problem(nep); }
741: if (jac) nep->computejacobian = jac;
742: if (ctx) nep->jacobianctx = ctx;
743: if (A) {
744: PetscObjectReference((PetscObject)A);
745: MatDestroy(&nep->jacobian);
746: nep->jacobian = A;
747: }
748: nep->fui = NEP_USER_INTERFACE_CALLBACK;
749: nep->state = NEP_STATE_INITIAL;
750: return(0);
751: }
753: /*@C
754: NEPGetJacobian - Returns the Jacobian matrix and optionally the user
755: provided routine and context for evaluating the Jacobian.
757: Not Collective, but Mat object will be parallel if NEP object is
759: Input Parameter:
760: . nep - the nonlinear eigensolver context
762: Output Parameters:
763: + A - location to stash Jacobian matrix (or NULL)
764: . jac - location to put Jacobian function (or NULL)
765: - ctx - location to stash Jacobian context (or NULL)
767: Level: advanced
769: .seealso: NEPSetJacobian()
770: @*/
771: PetscErrorCode NEPGetJacobian(NEP nep,Mat *A,PetscErrorCode (**jac)(NEP,PetscScalar,Mat,void*),void **ctx)772: {
775: NEPCheckCallback(nep,1);
776: if (A) *A = nep->jacobian;
777: if (jac) *jac = nep->computejacobian;
778: if (ctx) *ctx = nep->jacobianctx;
779: return(0);
780: }
782: /*@
783: NEPSetSplitOperator - Sets the operator of the nonlinear eigenvalue problem
784: in split form.
786: Collective on NEP, Mat and FN788: Input Parameters:
789: + nep - the nonlinear eigensolver context
790: . n - number of terms in the split form
791: . A - array of matrices
792: . f - array of functions
793: - str - structure flag for matrices
795: Notes:
796: The nonlinear operator is written as T(lambda) = sum_i A_i*f_i(lambda),
797: for i=1,...,n. The derivative T'(lambda) can be obtained using the
798: derivatives of f_i.
800: The structure flag provides information about A_i's nonzero pattern
801: (see MatStructure enum). If all matrices have the same pattern, then
802: use SAME_NONZERO_PATTERN. If the patterns are different but contained
803: in the pattern of the first one, then use SUBSET_NONZERO_PATTERN.
804: Otherwise use DIFFERENT_NONZERO_PATTERN.
806: This function must be called before NEPSetUp(). If it is called again
807: after NEPSetUp() then the NEP object is reset.
809: Level: beginner
811: .seealso: NEPGetSplitOperatorTerm(), NEPGetSplitOperatorInfo()
812: @*/
813: PetscErrorCode NEPSetSplitOperator(NEP nep,PetscInt n,Mat A[],FN f[],MatStructure str)814: {
815: PetscInt i;
821: if (n <= 0) SETERRQ1(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more terms, you have %D",n);
826: if (nep->state) { NEPReset(nep); }
827: else { NEPReset_Problem(nep); }
829: /* allocate space and copy matrices and functions */
830: PetscMalloc1(n,&nep->A);
831: PetscLogObjectMemory((PetscObject)nep,n*sizeof(Mat));
832: for (i=0;i<n;i++) {
834: PetscObjectReference((PetscObject)A[i]);
835: nep->A[i] = A[i];
836: }
837: PetscMalloc1(n,&nep->f);
838: PetscLogObjectMemory((PetscObject)nep,n*sizeof(FN));
839: for (i=0;i<n;i++) {
841: PetscObjectReference((PetscObject)f[i]);
842: nep->f[i] = f[i];
843: }
844: PetscCalloc1(n,&nep->nrma);
845: PetscLogObjectMemory((PetscObject)nep,n*sizeof(PetscReal));
846: nep->nt = n;
847: nep->mstr = str;
848: nep->fui = NEP_USER_INTERFACE_SPLIT;
849: nep->state = NEP_STATE_INITIAL;
850: return(0);
851: }
853: /*@
854: NEPGetSplitOperatorTerm - Gets the matrices and functions associated with
855: the nonlinear operator in split form.
857: Not collective, though parallel Mats and FNs are returned if the NEP is parallel
859: Input Parameter:
860: + nep - the nonlinear eigensolver context
861: - k - the index of the requested term (starting in 0)
863: Output Parameters:
864: + A - the matrix of the requested term
865: - f - the function of the requested term
867: Level: intermediate
869: .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorInfo()
870: @*/
871: PetscErrorCode NEPGetSplitOperatorTerm(NEP nep,PetscInt k,Mat *A,FN *f)872: {
875: NEPCheckSplit(nep,1);
876: if (k<0 || k>=nep->nt) SETERRQ1(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",nep->nt-1);
877: if (A) *A = nep->A[k];
878: if (f) *f = nep->f[k];
879: return(0);
880: }
882: /*@
883: NEPGetSplitOperatorInfo - Returns the number of terms of the split form of
884: the nonlinear operator, as well as the structure flag for matrices.
886: Not collective
888: Input Parameter:
889: . nep - the nonlinear eigensolver context
891: Output Parameters:
892: + n - the number of terms passed in NEPSetSplitOperator()
893: - str - the matrix structure flag passed in NEPSetSplitOperator()
895: Level: intermediate
897: .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorTerm()
898: @*/
899: PetscErrorCode NEPGetSplitOperatorInfo(NEP nep,PetscInt *n,MatStructure *str)900: {
903: NEPCheckSplit(nep,1);
904: if (n) *n = nep->nt;
905: if (str) *str = nep->mstr;
906: return(0);
907: }
909: /*@C
910: NEPSetDerivatives - Sets the function to compute the k-th derivative T^(k)(lambda)
911: for any value of k (including 0), as well as the location to store the matrix.
913: Logically Collective on NEP and Mat
915: Input Parameters:
916: + nep - the NEP context
917: . A - the matrix to store the computed derivative
918: . der - routing to evaluate the k-th derivative (if NULL then NEP retains any
919: previously set value)
920: - ctx - [optional] user-defined context for private data for the derivatives
921: evaluation routine (may be NULL) (if NULL then NEP retains any
922: previously set value)
924: Level: beginner
926: .seealso: NEPSetFunction(), NEPGetDerivatives()
927: @*/
928: PetscErrorCode NEPSetDerivatives(NEP nep,Mat A,PetscErrorCode (*der)(NEP,PetscScalar,PetscInt,Mat,void*),void *ctx)929: {
937: if (nep->state) { NEPReset(nep); }
938: else { NEPReset_Problem(nep); }
940: if (der) nep->computederivatives = der;
941: if (ctx) nep->derivativesctx = ctx;
942: if (A) {
943: PetscObjectReference((PetscObject)A);
944: MatDestroy(&nep->derivatives);
945: nep->derivatives = A;
946: }
947: nep->fui = NEP_USER_INTERFACE_DERIVATIVES;
948: nep->state = NEP_STATE_INITIAL;
949: return(0);
950: }
952: /*@C
953: NEPGetDerivatives - Returns the derivatives matrix and optionally the user
954: provided routine and context for evaluating the derivatives.
956: Not Collective, but Mat object will be parallel if NEP object is
958: Input Parameter:
959: . nep - the nonlinear eigensolver context
961: Output Parameters:
962: + A - location to stash the derivatives matrix (or NULL)
963: . der - location to put derivatives function (or NULL)
964: - ctx - location to stash derivatives context (or NULL)
966: Level: advanced
968: .seealso: NEPSetDerivatives()
969: @*/
970: PetscErrorCode NEPGetDerivatives(NEP nep,Mat *A,PetscErrorCode (**der)(NEP,PetscScalar,PetscInt,Mat,void*),void **ctx)971: {
974: NEPCheckDerivatives(nep,1);
975: if (A) *A = nep->derivatives;
976: if (der) *der = nep->computederivatives;
977: if (ctx) *ctx = nep->derivativesctx;
978: return(0);
979: }