1: /*
2: Basic NEP routines, Create, View, etc.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc-private/nepimpl.h> /*I "slepcnep.h" I*/
26: PetscFunctionList NEPList = 0;
27: PetscBool NEPRegisterAllCalled = PETSC_FALSE;
28: PetscClassId NEP_CLASSID = 0;
29: PetscLogEvent NEP_SetUp = 0,NEP_Solve = 0,NEP_Refine = 0,NEP_FunctionEval = 0,NEP_JacobianEval = 0;
33: /*@C
34: NEPView - Prints the NEP data structure.
36: Collective on NEP 38: Input Parameters:
39: + nep - the nonlinear eigenproblem solver context
40: - viewer - optional visualization context
42: Options Database Key:
43: . -nep_view - Calls NEPView() at end of NEPSolve()
45: Note:
46: The available visualization contexts include
47: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
48: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
49: output where only the first processor opens
50: the file. All other processors send their
51: data to the first processor to print.
53: The user can open an alternative visualization context with
54: PetscViewerASCIIOpen() - output to a specified file.
56: Level: beginner
58: .seealso: PetscViewerASCIIOpen()
59: @*/
60: PetscErrorCode NEPView(NEP nep,PetscViewer viewer) 61: {
63: char str[50];
64: PetscBool isascii,isslp,istrivial;
68: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)nep));
72: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
73: if (isascii) {
74: PetscObjectPrintClassNamePrefixType((PetscObject)nep,viewer);
75: if (nep->ops->view) {
76: PetscViewerASCIIPushTab(viewer);
77: (*nep->ops->view)(nep,viewer);
78: PetscViewerASCIIPopTab(viewer);
79: }
80: if (nep->split) {
81: PetscViewerASCIIPrintf(viewer," nonlinear operator in split form\n");
82: } else {
83: PetscViewerASCIIPrintf(viewer," nonlinear operator from user callbacks\n");
84: }
85: PetscViewerASCIIPrintf(viewer," iterative refinement: %s\n",NEPRefineTypes[nep->refine]);
86: if (nep->refine) {
87: PetscViewerASCIIPrintf(viewer," refinement stopping criterion: tol=%g, its=%D\n",(double)nep->reftol,nep->rits);
88: }
89: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: ");
90: SlepcSNPrintfScalar(str,50,nep->target,PETSC_FALSE);
91: if (!nep->which) {
92: PetscViewerASCIIPrintf(viewer,"not yet set\n");
93: } else switch (nep->which) {
94: case NEP_TARGET_MAGNITUDE:
95: PetscViewerASCIIPrintf(viewer,"closest to target: %s (in magnitude)\n",str);
96: break;
97: case NEP_TARGET_REAL:
98: PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the real axis)\n",str);
99: break;
100: case NEP_TARGET_IMAGINARY:
101: PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the imaginary axis)\n",str);
102: break;
103: case NEP_LARGEST_MAGNITUDE:
104: PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n");
105: break;
106: case NEP_SMALLEST_MAGNITUDE:
107: PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n");
108: break;
109: case NEP_LARGEST_REAL:
110: PetscViewerASCIIPrintf(viewer,"largest real parts\n");
111: break;
112: case NEP_SMALLEST_REAL:
113: PetscViewerASCIIPrintf(viewer,"smallest real parts\n");
114: break;
115: case NEP_LARGEST_IMAGINARY:
116: PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n");
117: break;
118: case NEP_SMALLEST_IMAGINARY:
119: PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n");
120: break;
121: default: SETERRQ(PetscObjectComm((PetscObject)nep),1,"Wrong value of nep->which");
122: }
123: PetscViewerASCIIPrintf(viewer," number of eigenvalues (nev): %D\n",nep->nev);
124: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %D\n",nep->ncv);
125: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %D\n",nep->mpd);
126: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %D\n",nep->max_it);
127: PetscViewerASCIIPrintf(viewer," maximum number of function evaluations: %D\n",nep->max_funcs);
128: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n",(double)nep->rtol,(double)nep->abstol,(double)nep->stol);
129: if (nep->lag) {
130: PetscViewerASCIIPrintf(viewer," updating the preconditioner every %D iterations\n",nep->lag);
131: }
132: if (nep->cctol) {
133: PetscViewerASCIIPrintf(viewer," using a constant tolerance for the linear solver\n");
134: }
135: if (nep->nini) {
136: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %D\n",PetscAbs(nep->nini));
137: }
138: } else {
139: if (nep->ops->view) {
140: (*nep->ops->view)(nep,viewer);
141: }
142: }
143: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
144: if (!nep->V) { NEPGetBV(nep,&nep->V); }
145: BVView(nep->V,viewer);
146: if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
147: RGIsTrivial(nep->rg,&istrivial);
148: if (!istrivial) { RGView(nep->rg,viewer); }
149: if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
150: DSView(nep->ds,viewer);
151: PetscViewerPopFormat(viewer);
152: PetscObjectTypeCompare((PetscObject)nep,NEPSLP,&isslp);
153: if (!isslp) {
154: if (!nep->ksp) { NEPGetKSP(nep,&nep->ksp); }
155: KSPView(nep->ksp,viewer);
156: }
157: return(0);
158: }
162: /*@C
163: NEPCreate - Creates the default NEP context.
165: Collective on MPI_Comm
167: Input Parameter:
168: . comm - MPI communicator
170: Output Parameter:
171: . nep - location to put the NEP context
173: Level: beginner
175: .seealso: NEPSetUp(), NEPSolve(), NEPDestroy(), NEP176: @*/
177: PetscErrorCode NEPCreate(MPI_Comm comm,NEP *outnep)178: {
180: NEP nep;
184: *outnep = 0;
185: NEPInitializePackage();
186: SlepcHeaderCreate(nep,_p_NEP,struct _NEPOps,NEP_CLASSID,"NEP","Nonlinear Eigenvalue Problem","NEP",comm,NEPDestroy,NEPView);
188: nep->max_it = 0;
189: nep->max_funcs = 0;
190: nep->nev = 1;
191: nep->ncv = 0;
192: nep->mpd = 0;
193: nep->lag = 1;
194: nep->nini = 0;
195: nep->target = 0.0;
196: nep->abstol = PETSC_DEFAULT;
197: nep->rtol = PETSC_DEFAULT;
198: nep->stol = PETSC_DEFAULT;
199: nep->ktol = 0.1;
200: nep->cctol = PETSC_FALSE;
201: nep->ttol = 0.0;
202: nep->which = (NEPWhich)0;
203: nep->refine = NEP_REFINE_NONE;
204: nep->reftol = PETSC_DEFAULT;
205: nep->rits = PETSC_DEFAULT;
206: nep->trackall = PETSC_FALSE;
208: nep->computefunction = NULL;
209: nep->computejacobian = NULL;
210: nep->functionctx = NULL;
211: nep->jacobianctx = NULL;
212: nep->converged = NEPConvergedDefault;
213: nep->convergeddestroy= NULL;
214: nep->convergedctx = NULL;
215: nep->numbermonitors = 0;
217: nep->ds = NULL;
218: nep->V = NULL;
219: nep->rg = NULL;
220: nep->rand = NULL;
221: nep->ksp = NULL;
222: nep->function = NULL;
223: nep->function_pre = NULL;
224: nep->jacobian = NULL;
225: nep->A = NULL;
226: nep->f = NULL;
227: nep->nt = 0;
228: nep->mstr = DIFFERENT_NONZERO_PATTERN;
229: nep->IS = NULL;
230: nep->eigr = NULL;
231: nep->eigi = NULL;
232: nep->errest = NULL;
233: nep->perm = NULL;
234: nep->nwork = 0;
235: nep->work = NULL;
236: nep->data = NULL;
238: nep->state = NEP_STATE_INITIAL;
239: nep->nconv = 0;
240: nep->its = 0;
241: nep->n = 0;
242: nep->nloc = 0;
243: nep->nfuncs = 0;
244: nep->split = PETSC_FALSE;
245: nep->reason = NEP_CONVERGED_ITERATING;
247: PetscNewLog(nep,&nep->sc);
248: PetscRandomCreate(comm,&nep->rand);
249: PetscRandomSetSeed(nep->rand,0x12345678);
250: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->rand);
251: *outnep = nep;
252: return(0);
253: }
257: /*@C
258: NEPSetType - Selects the particular solver to be used in the NEP object.
260: Logically Collective on NEP262: Input Parameters:
263: + nep - the nonlinear eigensolver context
264: - type - a known method
266: Options Database Key:
267: . -nep_type <method> - Sets the method; use -help for a list
268: of available methods
270: Notes:
271: See "slepc/include/slepcnep.h" for available methods.
273: Normally, it is best to use the NEPSetFromOptions() command and
274: then set the NEP type from the options database rather than by using
275: this routine. Using the options database provides the user with
276: maximum flexibility in evaluating the different available methods.
277: The NEPSetType() routine is provided for those situations where it
278: is necessary to set the iterative solver independently of the command
279: line or options database.
281: Level: intermediate
283: .seealso: NEPType284: @*/
285: PetscErrorCode NEPSetType(NEP nep,NEPType type)286: {
287: PetscErrorCode ierr,(*r)(NEP);
288: PetscBool match;
294: PetscObjectTypeCompare((PetscObject)nep,type,&match);
295: if (match) return(0);
297: PetscFunctionListFind(NEPList,type,&r);
298: if (!r) SETERRQ1(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown NEP type given: %s",type);
300: if (nep->ops->destroy) { (*nep->ops->destroy)(nep); }
301: PetscMemzero(nep->ops,sizeof(struct _NEPOps));
303: nep->state = NEP_STATE_INITIAL;
304: PetscObjectChangeTypeName((PetscObject)nep,type);
305: (*r)(nep);
306: return(0);
307: }
311: /*@C
312: NEPGetType - Gets the NEP type as a string from the NEP object.
314: Not Collective
316: Input Parameter:
317: . nep - the eigensolver context
319: Output Parameter:
320: . name - name of NEP method
322: Level: intermediate
324: .seealso: NEPSetType()
325: @*/
326: PetscErrorCode NEPGetType(NEP nep,NEPType *type)327: {
331: *type = ((PetscObject)nep)->type_name;
332: return(0);
333: }
337: /*@C
338: NEPRegister - Adds a method to the nonlinear eigenproblem solver package.
340: Not Collective
342: Input Parameters:
343: + name - name of a new user-defined solver
344: - function - routine to create the solver context
346: Notes:
347: NEPRegister() may be called multiple times to add several user-defined solvers.
349: Sample usage:
350: .vb
351: NEPRegister("my_solver",MySolverCreate);
352: .ve
354: Then, your solver can be chosen with the procedural interface via
355: $ NEPSetType(nep,"my_solver")
356: or at runtime via the option
357: $ -nep_type my_solver
359: Level: advanced
361: .seealso: NEPRegisterAll()
362: @*/
363: PetscErrorCode NEPRegister(const char *name,PetscErrorCode (*function)(NEP))364: {
368: PetscFunctionListAdd(&NEPList,name,function);
369: return(0);
370: }
374: /*@
375: NEPReset - Resets the NEP context to the initial state and removes any
376: allocated objects.
378: Collective on NEP380: Input Parameter:
381: . nep - eigensolver context obtained from NEPCreate()
383: Level: advanced
385: .seealso: NEPDestroy()
386: @*/
387: PetscErrorCode NEPReset(NEP nep)388: {
390: PetscInt i,ncols;
394: if (nep->ops->reset) { (nep->ops->reset)(nep); }
395: if (nep->ds) { DSReset(nep->ds); }
396: MatDestroy(&nep->function);
397: MatDestroy(&nep->function_pre);
398: MatDestroy(&nep->jacobian);
399: if (nep->split) {
400: MatDestroyMatrices(nep->nt,&nep->A);
401: for (i=0;i<nep->nt;i++) {
402: FNDestroy(&nep->f[i]);
403: }
404: PetscFree(nep->f);
405: }
406: BVGetSizes(nep->V,NULL,NULL,&ncols);
407: if (ncols) {
408: PetscFree4(nep->eigr,nep->eigi,nep->errest,nep->perm);
409: }
410: BVDestroy(&nep->V);
411: VecDestroyVecs(nep->nwork,&nep->work);
412: nep->nwork = 0;
413: nep->nfuncs = 0;
414: nep->state = NEP_STATE_INITIAL;
415: return(0);
416: }
420: /*@C
421: NEPDestroy - Destroys the NEP context.
423: Collective on NEP425: Input Parameter:
426: . nep - eigensolver context obtained from NEPCreate()
428: Level: beginner
430: .seealso: NEPCreate(), NEPSetUp(), NEPSolve()
431: @*/
432: PetscErrorCode NEPDestroy(NEP *nep)433: {
437: if (!*nep) return(0);
439: if (--((PetscObject)(*nep))->refct > 0) { *nep = 0; return(0); }
440: NEPReset(*nep);
441: if ((*nep)->ops->destroy) { (*(*nep)->ops->destroy)(*nep); }
442: KSPDestroy(&(*nep)->ksp);
443: RGDestroy(&(*nep)->rg);
444: DSDestroy(&(*nep)->ds);
445: PetscRandomDestroy(&(*nep)->rand);
446: PetscFree((*nep)->sc);
447: /* just in case the initial vectors have not been used */
448: SlepcBasisDestroy_Private(&(*nep)->nini,&(*nep)->IS);
449: if ((*nep)->convergeddestroy) {
450: (*(*nep)->convergeddestroy)((*nep)->convergedctx);
451: }
452: NEPMonitorCancel(*nep);
453: PetscHeaderDestroy(nep);
454: return(0);
455: }
459: /*@
460: NEPSetBV - Associates a basis vectors object to the nonlinear eigensolver.
462: Collective on NEP464: Input Parameters:
465: + nep - eigensolver context obtained from NEPCreate()
466: - bv - the basis vectors object
468: Note:
469: Use NEPGetBV() to retrieve the basis vectors context (for example,
470: to free it at the end of the computations).
472: Level: advanced
474: .seealso: NEPGetBV()
475: @*/
476: PetscErrorCode NEPSetBV(NEP nep,BV bv)477: {
484: PetscObjectReference((PetscObject)bv);
485: BVDestroy(&nep->V);
486: nep->V = bv;
487: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->V);
488: return(0);
489: }
493: /*@C
494: NEPGetBV - Obtain the basis vectors object associated to the nonlinear
495: eigensolver object.
497: Not Collective
499: Input Parameters:
500: . nep - eigensolver context obtained from NEPCreate()
502: Output Parameter:
503: . bv - basis vectors context
505: Level: advanced
507: .seealso: NEPSetBV()
508: @*/
509: PetscErrorCode NEPGetBV(NEP nep,BV *bv)510: {
516: if (!nep->V) {
517: BVCreate(PetscObjectComm((PetscObject)nep),&nep->V);
518: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->V);
519: }
520: *bv = nep->V;
521: return(0);
522: }
526: /*@
527: NEPSetRG - Associates a region object to the nonlinear eigensolver.
529: Collective on NEP531: Input Parameters:
532: + nep - eigensolver context obtained from NEPCreate()
533: - rg - the region object
535: Note:
536: Use NEPGetRG() to retrieve the region context (for example,
537: to free it at the end of the computations).
539: Level: advanced
541: .seealso: NEPGetRG()
542: @*/
543: PetscErrorCode NEPSetRG(NEP nep,RG rg)544: {
551: PetscObjectReference((PetscObject)rg);
552: RGDestroy(&nep->rg);
553: nep->rg = rg;
554: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->rg);
555: return(0);
556: }
560: /*@C
561: NEPGetRG - Obtain the region object associated to the
562: nonlinear eigensolver object.
564: Not Collective
566: Input Parameters:
567: . nep - eigensolver context obtained from NEPCreate()
569: Output Parameter:
570: . rg - region context
572: Level: advanced
574: .seealso: NEPSetRG()
575: @*/
576: PetscErrorCode NEPGetRG(NEP nep,RG *rg)577: {
583: if (!nep->rg) {
584: RGCreate(PetscObjectComm((PetscObject)nep),&nep->rg);
585: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->rg);
586: }
587: *rg = nep->rg;
588: return(0);
589: }
593: /*@
594: NEPSetDS - Associates a direct solver object to the nonlinear eigensolver.
596: Collective on NEP598: Input Parameters:
599: + nep - eigensolver context obtained from NEPCreate()
600: - ds - the direct solver object
602: Note:
603: Use NEPGetDS() to retrieve the direct solver context (for example,
604: to free it at the end of the computations).
606: Level: advanced
608: .seealso: NEPGetDS()
609: @*/
610: PetscErrorCode NEPSetDS(NEP nep,DS ds)611: {
618: PetscObjectReference((PetscObject)ds);
619: DSDestroy(&nep->ds);
620: nep->ds = ds;
621: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->ds);
622: return(0);
623: }
627: /*@C
628: NEPGetDS - Obtain the direct solver object associated to the
629: nonlinear eigensolver object.
631: Not Collective
633: Input Parameters:
634: . nep - eigensolver context obtained from NEPCreate()
636: Output Parameter:
637: . ds - direct solver context
639: Level: advanced
641: .seealso: NEPSetDS()
642: @*/
643: PetscErrorCode NEPGetDS(NEP nep,DS *ds)644: {
650: if (!nep->ds) {
651: DSCreate(PetscObjectComm((PetscObject)nep),&nep->ds);
652: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->ds);
653: }
654: *ds = nep->ds;
655: return(0);
656: }
660: /*@
661: NEPSetKSP - Associates a linear solver object to the nonlinear eigensolver.
663: Collective on NEP665: Input Parameters:
666: + nep - eigensolver context obtained from NEPCreate()
667: - ksp - the linear solver object
669: Note:
670: Use NEPGetKSP() to retrieve the linear solver context (for example,
671: to free it at the end of the computations).
673: Level: developer
675: .seealso: NEPGetKSP()
676: @*/
677: PetscErrorCode NEPSetKSP(NEP nep,KSP ksp)678: {
685: PetscObjectReference((PetscObject)ksp);
686: KSPDestroy(&nep->ksp);
687: nep->ksp = ksp;
688: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->ksp);
689: return(0);
690: }
694: /*@C
695: NEPGetKSP - Obtain the linear solver (KSP) object associated
696: to the eigensolver object.
698: Not Collective
700: Input Parameters:
701: . nep - eigensolver context obtained from NEPCreate()
703: Output Parameter:
704: . ksp - linear solver context
706: Level: beginner
708: .seealso: NEPSetKSP()
709: @*/
710: PetscErrorCode NEPGetKSP(NEP nep,KSP *ksp)711: {
717: if (!nep->ksp) {
718: KSPCreate(PetscObjectComm((PetscObject)nep),&nep->ksp);
719: KSPSetOptionsPrefix(nep->ksp,((PetscObject)nep)->prefix);
720: KSPAppendOptionsPrefix(nep->ksp,"nep_");
721: PetscObjectIncrementTabLevel((PetscObject)nep->ksp,(PetscObject)nep,1);
722: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->ksp);
723: }
724: *ksp = nep->ksp;
725: return(0);
726: }
730: /*@
731: NEPSetTarget - Sets the value of the target.
733: Logically Collective on NEP735: Input Parameters:
736: + nep - eigensolver context
737: - target - the value of the target
739: Options Database Key:
740: . -nep_target <scalar> - the value of the target
742: Notes:
743: The target is a scalar value used to determine the portion of the spectrum
744: of interest. It is used in combination with NEPSetWhichEigenpairs().
746: In the case of complex scalars, a complex value can be provided in the
747: command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
748: -nep_target 1.0+2.0i
750: Level: beginner
752: .seealso: NEPGetTarget(), NEPSetWhichEigenpairs()
753: @*/
754: PetscErrorCode NEPSetTarget(NEP nep,PetscScalar target)755: {
759: nep->target = target;
760: return(0);
761: }
765: /*@
766: NEPGetTarget - Gets the value of the target.
768: Not Collective
770: Input Parameter:
771: . nep - eigensolver context
773: Output Parameter:
774: . target - the value of the target
776: Note:
777: If the target was not set by the user, then zero is returned.
779: Level: beginner
781: .seealso: NEPSetTarget()
782: @*/
783: PetscErrorCode NEPGetTarget(NEP nep,PetscScalar* target)784: {
788: *target = nep->target;
789: return(0);
790: }
794: /*@C
795: NEPSetFunction - Sets the function to compute the nonlinear Function T(lambda)
796: as well as the location to store the matrix.
798: Logically Collective on NEP and Mat
800: Input Parameters:
801: + nep - the NEP context
802: . A - Function matrix
803: . B - preconditioner matrix (usually same as the Function)
804: . fun - Function evaluation routine (if NULL then NEP retains any
805: previously set value)
806: - ctx - [optional] user-defined context for private data for the Function
807: evaluation routine (may be NULL) (if NULL then NEP retains any
808: previously set value)
810: Level: beginner
812: .seealso: NEPGetFunction(), NEPSetJacobian()
813: @*/
814: PetscErrorCode NEPSetFunction(NEP nep,Mat A,Mat B,PetscErrorCode (*fun)(NEP,PetscScalar,Mat,Mat,void*),void *ctx)815: {
824: if (fun) nep->computefunction = fun;
825: if (ctx) nep->functionctx = ctx;
826: if (A) {
827: PetscObjectReference((PetscObject)A);
828: MatDestroy(&nep->function);
829: nep->function = A;
830: }
831: if (B) {
832: PetscObjectReference((PetscObject)B);
833: MatDestroy(&nep->function_pre);
834: nep->function_pre = B;
835: }
836: nep->split = PETSC_FALSE;
837: return(0);
838: }
842: /*@C
843: NEPGetFunction - Returns the Function matrix and optionally the user
844: provided context for evaluating the Function.
846: Not Collective, but Mat object will be parallel if NEP object is
848: Input Parameter:
849: . nep - the nonlinear eigensolver context
851: Output Parameters:
852: + A - location to stash Function matrix (or NULL)
853: . B - location to stash preconditioner matrix (or NULL)
854: . fun - location to put Function function (or NULL)
855: - ctx - location to stash Function context (or NULL)
857: Level: advanced
859: .seealso: NEPSetFunction()
860: @*/
861: PetscErrorCode NEPGetFunction(NEP nep,Mat *A,Mat *B,PetscErrorCode (**fun)(NEP,PetscScalar,Mat,Mat,void*),void **ctx)862: {
865: if (A) *A = nep->function;
866: if (B) *B = nep->function_pre;
867: if (fun) *fun = nep->computefunction;
868: if (ctx) *ctx = nep->functionctx;
869: return(0);
870: }
874: /*@C
875: NEPSetJacobian - Sets the function to compute Jacobian T'(lambda) as well
876: as the location to store the matrix.
878: Logically Collective on NEP and Mat
880: Input Parameters:
881: + nep - the NEP context
882: . A - Jacobian matrix
883: . jac - Jacobian evaluation routine (if NULL then NEP retains any
884: previously set value)
885: - ctx - [optional] user-defined context for private data for the Jacobian
886: evaluation routine (may be NULL) (if NULL then NEP retains any
887: previously set value)
889: Level: beginner
891: .seealso: NEPSetFunction(), NEPGetJacobian()
892: @*/
893: PetscErrorCode NEPSetJacobian(NEP nep,Mat A,PetscErrorCode (*jac)(NEP,PetscScalar,Mat,void*),void *ctx)894: {
901: if (jac) nep->computejacobian = jac;
902: if (ctx) nep->jacobianctx = ctx;
903: if (A) {
904: PetscObjectReference((PetscObject)A);
905: MatDestroy(&nep->jacobian);
906: nep->jacobian = A;
907: }
908: nep->split = PETSC_FALSE;
909: return(0);
910: }
914: /*@C
915: NEPGetJacobian - Returns the Jacobian matrix and optionally the user
916: provided context for evaluating the Jacobian.
918: Not Collective, but Mat object will be parallel if NEP object is
920: Input Parameter:
921: . nep - the nonlinear eigensolver context
923: Output Parameters:
924: + A - location to stash Jacobian matrix (or NULL)
925: . jac - location to put Jacobian function (or NULL)
926: - ctx - location to stash Jacobian context (or NULL)
928: Level: advanced
930: .seealso: NEPSetJacobian()
931: @*/
932: PetscErrorCode NEPGetJacobian(NEP nep,Mat *A,PetscErrorCode (**jac)(NEP,PetscScalar,Mat,void*),void **ctx)933: {
936: if (A) *A = nep->jacobian;
937: if (jac) *jac = nep->computejacobian;
938: if (ctx) *ctx = nep->jacobianctx;
939: return(0);
940: }
944: /*@
945: NEPSetSplitOperator - Sets the operator of the nonlinear eigenvalue problem
946: in split form.
948: Collective on NEP, Mat and FN950: Input Parameters:
951: + nep - the nonlinear eigensolver context
952: . n - number of terms in the split form
953: . A - array of matrices
954: . f - array of functions
955: - str - structure flag for matrices
957: Notes:
958: The nonlinear operator is written as T(lambda) = sum_i A_i*f_i(lambda),
959: for i=1,...,n. The derivative T'(lambda) can be obtained using the
960: derivatives of f_i.
962: The structure flag provides information about A_i's nonzero pattern
963: (see MatStructure enum). If all matrices have the same pattern, then
964: use SAME_NONZERO_PATTERN. If the patterns are different but contained
965: in the pattern of the first one, then use SUBSET_NONZERO_PATTERN.
966: Otherwise use DIFFERENT_NONZERO_PATTERN.
968: This function must be called before NEPSetUp(). If it is called again
969: after NEPSetUp() then the NEP object is reset.
971: Level: intermediate
973: .seealso: NEPGetSplitOperatorTerm(), NEPGetSplitOperatorInfo()
974: @*/
975: PetscErrorCode NEPSetSplitOperator(NEP nep,PetscInt n,Mat A[],FN f[],MatStructure str)976: {
977: PetscInt i;
983: if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more terms, you have %D",n);
988: if (nep->state) { NEPReset(nep); }
989: /* clean previously stored information */
990: MatDestroy(&nep->function);
991: MatDestroy(&nep->function_pre);
992: MatDestroy(&nep->jacobian);
993: if (nep->split) {
994: MatDestroyMatrices(nep->nt,&nep->A);
995: for (i=0;i<nep->nt;i++) {
996: FNDestroy(&nep->f[i]);
997: }
998: PetscFree(nep->f);
999: }
1000: /* allocate space and copy matrices and functions */
1001: PetscMalloc1(n,&nep->A);
1002: PetscLogObjectMemory((PetscObject)nep,n*sizeof(Mat));
1003: for (i=0;i<n;i++) {
1005: PetscObjectReference((PetscObject)A[i]);
1006: nep->A[i] = A[i];
1007: }
1008: PetscMalloc1(n,&nep->f);
1009: PetscLogObjectMemory((PetscObject)nep,n*sizeof(FN));
1010: for (i=0;i<n;i++) {
1012: PetscObjectReference((PetscObject)f[i]);
1013: nep->f[i] = f[i];
1014: }
1015: nep->nt = n;
1016: nep->mstr = str;
1017: nep->split = PETSC_TRUE;
1018: return(0);
1019: }
1023: /*@
1024: NEPGetSplitOperatorTerm - Gets the matrices and functions associated with
1025: the nonlinear operator in split form.
1027: Not collective, though parallel Mats and FNs are returned if the NEP is parallel
1029: Input Parameter:
1030: + nep - the nonlinear eigensolver context
1031: - k - the index of the requested term (starting in 0)
1033: Output Parameters:
1034: + A - the matrix of the requested term
1035: - f - the function of the requested term
1037: Level: intermediate
1039: .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorInfo()
1040: @*/
1041: PetscErrorCode NEPGetSplitOperatorTerm(NEP nep,PetscInt k,Mat *A,FN *f)1042: {
1045: if (k<0 || k>=nep->nt) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %d",nep->nt-1);
1046: if (A) *A = nep->A[k];
1047: if (f) *f = nep->f[k];
1048: return(0);
1049: }
1053: /*@
1054: NEPGetSplitOperatorInfo - Returns the number of terms of the split form of
1055: the nonlinear operator, as well as the structure flag for matrices.
1057: Not collective
1059: Input Parameter:
1060: . nep - the nonlinear eigensolver context
1062: Output Parameters:
1063: + n - the number of terms passed in NEPSetSplitOperator()
1064: - str - the matrix structure flag passed in NEPSetSplitOperator()
1066: Level: intermediate
1068: .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorTerm()
1069: @*/
1070: PetscErrorCode NEPGetSplitOperatorInfo(NEP nep,PetscInt *n,MatStructure *str)1071: {
1074: if (n) *n = nep->nt;
1075: if (str) *str = nep->mstr;
1076: return(0);
1077: }