Actual source code: nepbasic.c
1: /*
2: Basic NEP routines, Create, View, etc.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2013, 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_Dense = 0,NEP_FunctionEval = 0,NEP_JacobianEval = 0;
30: static PetscBool NEPPackageInitialized = PETSC_FALSE;
34: /*@C
35: NEPFinalizePackage - This function destroys everything in the Slepc interface
36: to the NEP package. It is called from SlepcFinalize().
38: Level: developer
40: .seealso: SlepcFinalize()
41: @*/
42: PetscErrorCode NEPFinalizePackage(void)
43: {
47: PetscFunctionListDestroy(&NEPList);
48: NEPPackageInitialized = PETSC_FALSE;
49: NEPRegisterAllCalled = PETSC_FALSE;
50: return(0);
51: }
55: /*@C
56: NEPInitializePackage - This function initializes everything in the NEP package. It is called
57: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to NEPCreate()
58: when using static libraries.
60: Level: developer
62: .seealso: SlepcInitialize()
63: @*/
64: PetscErrorCode NEPInitializePackage(void)
65: {
66: char logList[256];
67: char *className;
68: PetscBool opt;
72: if (NEPPackageInitialized) return(0);
73: NEPPackageInitialized = PETSC_TRUE;
74: /* Register Classes */
75: PetscClassIdRegister("Nonlinear Eigenvalue Problem solver",&NEP_CLASSID);
76: /* Register Constructors */
77: NEPRegisterAll();
78: /* Register Events */
79: PetscLogEventRegister("NEPSetUp",NEP_CLASSID,&NEP_SetUp);
80: PetscLogEventRegister("NEPSolve",NEP_CLASSID,&NEP_Solve);
81: PetscLogEventRegister("NEPDense",NEP_CLASSID,&NEP_Dense);
82: PetscLogEventRegister("NEPFunctionEval",NEP_CLASSID,&NEP_FunctionEval);
83: PetscLogEventRegister("NEPJacobianEval",NEP_CLASSID,&NEP_JacobianEval);
84: /* Process info exclusions */
85: PetscOptionsGetString(NULL,"-info_exclude",logList,256,&opt);
86: if (opt) {
87: PetscStrstr(logList,"nep",&className);
88: if (className) {
89: PetscInfoDeactivateClass(NEP_CLASSID);
90: }
91: }
92: /* Process summary exclusions */
93: PetscOptionsGetString(NULL,"-log_summary_exclude",logList,256,&opt);
94: if (opt) {
95: PetscStrstr(logList,"nep",&className);
96: if (className) {
97: PetscLogEventDeactivateClass(NEP_CLASSID);
98: }
99: }
100: PetscRegisterFinalize(NEPFinalizePackage);
101: return(0);
102: }
106: /*@C
107: NEPView - Prints the NEP data structure.
109: Collective on NEP
111: Input Parameters:
112: + nep - the nonlinear eigenproblem solver context
113: - viewer - optional visualization context
115: Options Database Key:
116: . -nep_view - Calls NEPView() at end of NEPSolve()
118: Note:
119: The available visualization contexts include
120: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
121: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
122: output where only the first processor opens
123: the file. All other processors send their
124: data to the first processor to print.
126: The user can open an alternative visualization context with
127: PetscViewerASCIIOpen() - output to a specified file.
129: Level: beginner
131: .seealso: PetscViewerASCIIOpen()
132: @*/
133: PetscErrorCode NEPView(NEP nep,PetscViewer viewer)
134: {
136: char str[50];
137: PetscBool isascii,isslp;
141: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)nep));
145: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
146: if (isascii) {
147: PetscObjectPrintClassNamePrefixType((PetscObject)nep,viewer,"NEP Object");
148: if (nep->ops->view) {
149: PetscViewerASCIIPushTab(viewer);
150: (*nep->ops->view)(nep,viewer);
151: PetscViewerASCIIPopTab(viewer);
152: }
153: if (nep->split) {
154: PetscViewerASCIIPrintf(viewer," nonlinear operator in split form\n");
155: } else {
156: PetscViewerASCIIPrintf(viewer," nonlinear operator from user callbacks\n");
157: }
158: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: ");
159: SlepcSNPrintfScalar(str,50,nep->target,PETSC_FALSE);
160: if (!nep->which) {
161: PetscViewerASCIIPrintf(viewer,"not yet set\n");
162: } else switch (nep->which) {
163: case NEP_TARGET_MAGNITUDE:
164: PetscViewerASCIIPrintf(viewer,"closest to target: %s (in magnitude)\n",str);
165: break;
166: case NEP_TARGET_REAL:
167: PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the real axis)\n",str);
168: break;
169: case NEP_TARGET_IMAGINARY:
170: PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the imaginary axis)\n",str);
171: break;
172: case NEP_LARGEST_MAGNITUDE:
173: PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n");
174: break;
175: case NEP_SMALLEST_MAGNITUDE:
176: PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n");
177: break;
178: case NEP_LARGEST_REAL:
179: PetscViewerASCIIPrintf(viewer,"largest real parts\n");
180: break;
181: case NEP_SMALLEST_REAL:
182: PetscViewerASCIIPrintf(viewer,"smallest real parts\n");
183: break;
184: case NEP_LARGEST_IMAGINARY:
185: PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n");
186: break;
187: case NEP_SMALLEST_IMAGINARY:
188: PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n");
189: break;
190: default: SETERRQ(PetscObjectComm((PetscObject)nep),1,"Wrong value of nep->which");
191: }
192: PetscViewerASCIIPrintf(viewer," number of eigenvalues (nev): %D\n",nep->nev);
193: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %D\n",nep->ncv);
194: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %D\n",nep->mpd);
195: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %D\n",nep->max_it);
196: PetscViewerASCIIPrintf(viewer," maximum number of function evaluations: %D\n",nep->max_funcs);
197: PetscViewerASCIIPrintf(viewer," tolerances: relative=%G, absolute=%G, solution=%G\n",nep->rtol,nep->abstol,nep->stol);
198: if (nep->lag) {
199: PetscViewerASCIIPrintf(viewer," updating the preconditioner every %D iterations\n",nep->lag);
200: }
201: if (nep->cctol) {
202: PetscViewerASCIIPrintf(viewer," using a constant tolerance for the linear solver\n");
203: }
204: if (nep->nini) {
205: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %D\n",PetscAbs(nep->nini));
206: }
207: } else {
208: if (nep->ops->view) {
209: (*nep->ops->view)(nep,viewer);
210: }
211: }
212: if (!nep->ip) { NEPGetIP(nep,&nep->ip); }
213: IPView(nep->ip,viewer);
214: if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
215: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
216: DSView(nep->ds,viewer);
217: PetscViewerPopFormat(viewer);
218: PetscObjectTypeCompare((PetscObject)nep,NEPSLP,&isslp);
219: if (!isslp) {
220: if (!nep->ksp) { NEPGetKSP(nep,&nep->ksp); }
221: KSPView(nep->ksp,viewer);
222: }
223: return(0);
224: }
228: /*@C
229: NEPCreate - Creates the default NEP context.
231: Collective on MPI_Comm
233: Input Parameter:
234: . comm - MPI communicator
236: Output Parameter:
237: . nep - location to put the NEP context
239: Level: beginner
241: .seealso: NEPSetUp(), NEPSolve(), NEPDestroy(), NEP
242: @*/
243: PetscErrorCode NEPCreate(MPI_Comm comm,NEP *outnep)
244: {
246: NEP nep;
250: *outnep = 0;
251: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
252: NEPInitializePackage();
253: #endif
255: SlepcHeaderCreate(nep,_p_NEP,struct _NEPOps,NEP_CLASSID,"NEP","Nonlinear Eigenvalue Problem","NEP",comm,NEPDestroy,NEPView);
257: nep->max_it = 0;
258: nep->max_funcs = 0;
259: nep->nev = 1;
260: nep->ncv = 0;
261: nep->mpd = 0;
262: nep->lag = 1;
263: nep->nini = 0;
264: nep->allocated_ncv = 0;
265: nep->ip = 0;
266: nep->ds = 0;
267: nep->function = 0;
268: nep->function_pre = 0;
269: nep->jacobian = 0;
270: nep->abstol = PETSC_DEFAULT;
271: nep->rtol = PETSC_DEFAULT;
272: nep->stol = PETSC_DEFAULT;
273: nep->ktol = 0.1;
274: nep->cctol = PETSC_FALSE;
275: nep->ttol = 0.0;
276: nep->which = (NEPWhich)0;
277: nep->computefunction = NULL;
278: nep->computejacobian = NULL;
279: nep->comparison = NULL;
280: nep->converged = NEPConvergedDefault;
281: nep->convergeddestroy= NULL;
282: nep->comparisonctx = NULL;
283: nep->convergedctx = NULL;
284: nep->functionctx = NULL;
285: nep->jacobianctx = NULL;
286: nep->V = NULL;
287: nep->IS = NULL;
288: nep->eig = NULL;
289: nep->errest = NULL;
290: nep->data = NULL;
291: nep->t = NULL;
292: nep->split = PETSC_FALSE;
293: nep->nt = 0;
294: nep->mstr = DIFFERENT_NONZERO_PATTERN;
295: nep->A = NULL;
296: nep->f = NULL;
297: nep->nconv = 0;
298: nep->its = 0;
299: nep->perm = NULL;
300: nep->nfuncs = 0;
301: nep->linits = 0;
302: nep->nwork = 0;
303: nep->work = NULL;
304: nep->setupcalled = 0;
305: nep->reason = NEP_CONVERGED_ITERATING;
306: nep->numbermonitors = 0;
307: nep->trackall = PETSC_FALSE;
308: nep->rand = 0;
310: PetscRandomCreate(comm,&nep->rand);
311: PetscRandomSetSeed(nep->rand,0x12345678);
312: PetscLogObjectParent(nep,nep->rand);
313: *outnep = nep;
314: return(0);
315: }
319: /*@C
320: NEPSetType - Selects the particular solver to be used in the NEP object.
322: Logically Collective on NEP
324: Input Parameters:
325: + nep - the nonlinear eigensolver context
326: - type - a known method
328: Options Database Key:
329: . -nep_type <method> - Sets the method; use -help for a list
330: of available methods
332: Notes:
333: See "slepc/include/slepcnep.h" for available methods.
335: Normally, it is best to use the NEPSetFromOptions() command and
336: then set the NEP type from the options database rather than by using
337: this routine. Using the options database provides the user with
338: maximum flexibility in evaluating the different available methods.
339: The NEPSetType() routine is provided for those situations where it
340: is necessary to set the iterative solver independently of the command
341: line or options database.
343: Level: intermediate
345: .seealso: NEPType
346: @*/
347: PetscErrorCode NEPSetType(NEP nep,NEPType type)
348: {
349: PetscErrorCode ierr,(*r)(NEP);
350: PetscBool match;
356: PetscObjectTypeCompare((PetscObject)nep,type,&match);
357: if (match) return(0);
359: PetscFunctionListFind(NEPList,type,&r);
360: if (!r) SETERRQ1(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown NEP type given: %s",type);
362: if (nep->ops->destroy) { (*nep->ops->destroy)(nep); }
363: PetscMemzero(nep->ops,sizeof(struct _NEPOps));
365: nep->setupcalled = 0;
366: PetscObjectChangeTypeName((PetscObject)nep,type);
367: (*r)(nep);
368: return(0);
369: }
373: /*@C
374: NEPGetType - Gets the NEP type as a string from the NEP object.
376: Not Collective
378: Input Parameter:
379: . nep - the eigensolver context
381: Output Parameter:
382: . name - name of NEP method
384: Level: intermediate
386: .seealso: NEPSetType()
387: @*/
388: PetscErrorCode NEPGetType(NEP nep,NEPType *type)
389: {
393: *type = ((PetscObject)nep)->type_name;
394: return(0);
395: }
399: /*@C
400: NEPRegister - Adds a method to the quadratic eigenproblem solver package.
402: Not Collective
404: Input Parameters:
405: + name - name of a new user-defined solver
406: - function - routine to create the solver context
408: Notes:
409: NEPRegister() may be called multiple times to add several user-defined solvers.
411: Sample usage:
412: .vb
413: NEPRegister("my_solver",MySolverCreate);
414: .ve
416: Then, your solver can be chosen with the procedural interface via
417: $ NEPSetType(qep,"my_solver")
418: or at runtime via the option
419: $ -qep_type my_solver
421: Level: advanced
423: .seealso: NEPRegisterAll()
424: @*/
425: PetscErrorCode NEPRegister(const char *name,PetscErrorCode (*function)(NEP))
426: {
430: PetscFunctionListAdd(&NEPList,name,function);
431: return(0);
432: }
436: /*@
437: NEPReset - Resets the NEP context to the setupcalled=0 state and removes any
438: allocated objects.
440: Collective on NEP
442: Input Parameter:
443: . nep - eigensolver context obtained from NEPCreate()
445: Level: advanced
447: .seealso: NEPDestroy()
448: @*/
449: PetscErrorCode NEPReset(NEP nep)
450: {
455: if (nep->ops->reset) { (nep->ops->reset)(nep); }
456: if (nep->ip) { IPReset(nep->ip); }
457: if (nep->ds) { DSReset(nep->ds); }
458: VecDestroy(&nep->t);
459: NEPFreeSolution(nep);
460: nep->nfuncs = 0;
461: nep->linits = 0;
462: nep->setupcalled = 0;
463: return(0);
464: }
468: /*@C
469: NEPDestroy - Destroys the NEP context.
471: Collective on NEP
473: Input Parameter:
474: . nep - eigensolver context obtained from NEPCreate()
476: Level: beginner
478: .seealso: NEPCreate(), NEPSetUp(), NEPSolve()
479: @*/
480: PetscErrorCode NEPDestroy(NEP *nep)
481: {
483: PetscInt i;
486: if (!*nep) return(0);
488: if (--((PetscObject)(*nep))->refct > 0) { *nep = 0; return(0); }
489: NEPReset(*nep);
490: if ((*nep)->ops->destroy) { (*(*nep)->ops->destroy)(*nep); }
491: KSPDestroy(&(*nep)->ksp);
492: IPDestroy(&(*nep)->ip);
493: DSDestroy(&(*nep)->ds);
494: MatDestroy(&(*nep)->function);
495: MatDestroy(&(*nep)->function_pre);
496: MatDestroy(&(*nep)->jacobian);
497: if ((*nep)->split) {
498: MatDestroyMatrices((*nep)->nt,&(*nep)->A);
499: for (i=0;i<(*nep)->nt;i++) {
500: FNDestroy(&(*nep)->f[i]);
501: }
502: PetscFree((*nep)->f);
503: }
504: PetscRandomDestroy(&(*nep)->rand);
505: /* just in case the initial vectors have not been used */
506: SlepcBasisDestroy_Private(&(*nep)->nini,&(*nep)->IS);
507: NEPMonitorCancel(*nep);
508: PetscHeaderDestroy(nep);
509: return(0);
510: }
514: /*@
515: NEPSetIP - Associates an inner product object to the nonlinear eigensolver.
517: Collective on NEP
519: Input Parameters:
520: + nep - eigensolver context obtained from NEPCreate()
521: - ip - the inner product object
523: Note:
524: Use NEPGetIP() to retrieve the inner product context (for example,
525: to free it at the end of the computations).
527: Level: advanced
529: .seealso: NEPGetIP()
530: @*/
531: PetscErrorCode NEPSetIP(NEP nep,IP ip)
532: {
539: PetscObjectReference((PetscObject)ip);
540: IPDestroy(&nep->ip);
541: nep->ip = ip;
542: PetscLogObjectParent(nep,nep->ip);
543: return(0);
544: }
548: /*@C
549: NEPGetIP - Obtain the inner product object associated
550: to the nonlinear eigensolver object.
552: Not Collective
554: Input Parameters:
555: . nep - eigensolver context obtained from NEPCreate()
557: Output Parameter:
558: . ip - inner product context
560: Level: advanced
562: .seealso: NEPSetIP()
563: @*/
564: PetscErrorCode NEPGetIP(NEP nep,IP *ip)
565: {
571: if (!nep->ip) {
572: IPCreate(PetscObjectComm((PetscObject)nep),&nep->ip);
573: PetscLogObjectParent(nep,nep->ip);
574: }
575: *ip = nep->ip;
576: return(0);
577: }
581: /*@
582: NEPSetDS - Associates a direct solver object to the nonlinear eigensolver.
584: Collective on NEP
586: Input Parameters:
587: + nep - eigensolver context obtained from NEPCreate()
588: - ds - the direct solver object
590: Note:
591: Use NEPGetDS() to retrieve the direct solver context (for example,
592: to free it at the end of the computations).
594: Level: advanced
596: .seealso: NEPGetDS()
597: @*/
598: PetscErrorCode NEPSetDS(NEP nep,DS ds)
599: {
606: PetscObjectReference((PetscObject)ds);
607: DSDestroy(&nep->ds);
608: nep->ds = ds;
609: PetscLogObjectParent(nep,nep->ds);
610: return(0);
611: }
615: /*@C
616: NEPGetDS - Obtain the direct solver object associated to the
617: nonlinear eigensolver object.
619: Not Collective
621: Input Parameters:
622: . nep - eigensolver context obtained from NEPCreate()
624: Output Parameter:
625: . ds - direct solver context
627: Level: advanced
629: .seealso: NEPSetDS()
630: @*/
631: PetscErrorCode NEPGetDS(NEP nep,DS *ds)
632: {
638: if (!nep->ds) {
639: DSCreate(PetscObjectComm((PetscObject)nep),&nep->ds);
640: PetscLogObjectParent(nep,nep->ds);
641: }
642: *ds = nep->ds;
643: return(0);
644: }
648: /*@
649: NEPSetKSP - Associates a linear solver object to the nonlinear eigensolver.
651: Collective on NEP
653: Input Parameters:
654: + nep - eigensolver context obtained from NEPCreate()
655: - ksp - the linear solver object
657: Note:
658: Use NEPGetKSP() to retrieve the linear solver context (for example,
659: to free it at the end of the computations).
661: Level: developer
663: .seealso: NEPGetKSP()
664: @*/
665: PetscErrorCode NEPSetKSP(NEP nep,KSP ksp)
666: {
673: PetscObjectReference((PetscObject)ksp);
674: KSPDestroy(&nep->ksp);
675: nep->ksp = ksp;
676: PetscLogObjectParent(nep,nep->ksp);
677: return(0);
678: }
682: /*@C
683: NEPGetKSP - Obtain the linear solver (KSP) object associated
684: to the eigensolver object.
686: Not Collective
688: Input Parameters:
689: . nep - eigensolver context obtained from NEPCreate()
691: Output Parameter:
692: . ksp - linear solver context
694: Level: beginner
696: .seealso: NEPSetKSP()
697: @*/
698: PetscErrorCode NEPGetKSP(NEP nep,KSP *ksp)
699: {
705: if (!nep->ksp) {
706: KSPCreate(PetscObjectComm((PetscObject)nep),&nep->ksp);
707: KSPSetOptionsPrefix(nep->ksp,((PetscObject)nep)->prefix);
708: KSPAppendOptionsPrefix(nep->ksp,"nep_");
709: PetscObjectIncrementTabLevel((PetscObject)nep->ksp,(PetscObject)nep,1);
710: PetscLogObjectParent(nep,nep->ksp);
711: }
712: *ksp = nep->ksp;
713: return(0);
714: }
718: /*@
719: NEPSetTarget - Sets the value of the target.
721: Logically Collective on NEP
723: Input Parameters:
724: + nep - eigensolver context
725: - target - the value of the target
727: Notes:
728: The target is a scalar value used to determine the portion of the spectrum
729: of interest. It is used in combination with NEPSetWhichEigenpairs().
731: Level: beginner
733: .seealso: NEPGetTarget(), NEPSetWhichEigenpairs()
734: @*/
735: PetscErrorCode NEPSetTarget(NEP nep,PetscScalar target)
736: {
740: nep->target = target;
741: return(0);
742: }
746: /*@
747: NEPGetTarget - Gets the value of the target.
749: Not Collective
751: Input Parameter:
752: . nep - eigensolver context
754: Output Parameter:
755: . target - the value of the target
757: Level: beginner
759: Note:
760: If the target was not set by the user, then zero is returned.
762: .seealso: NEPSetTarget()
763: @*/
764: PetscErrorCode NEPGetTarget(NEP nep,PetscScalar* target)
765: {
769: *target = nep->target;
770: return(0);
771: }
775: /*@C
776: NEPSetFunction - Sets the function to compute the nonlinear Function T(lambda)
777: as well as the location to store the matrix.
779: Logically Collective on NEP and Mat
781: Input Parameters:
782: + nep - the NEP context
783: . A - Function matrix
784: . B - preconditioner matrix (usually same as the Function)
785: . fun - Function evaluation routine (if NULL then NEP retains any
786: previously set value)
787: - ctx - [optional] user-defined context for private data for the Function
788: evaluation routine (may be NULL) (if NULL then NEP retains any
789: previously set value)
791: Notes:
792: The routine fun() takes Mat* as the matrix arguments rather than Mat.
793: This allows the Function evaluation routine to replace A and/or B with a
794: completely new matrix structure (not just different matrix elements)
795: when appropriate, for instance, if the nonzero structure is changing
796: throughout the global iterations.
798: Level: beginner
800: .seealso: NEPGetFunction(), NEPSetJacobian()
801: @*/
802: PetscErrorCode NEPSetFunction(NEP nep,Mat A,Mat B,PetscErrorCode (*fun)(NEP,PetscScalar,Mat*,Mat*,MatStructure*,void*),void *ctx)
803: {
812: if (fun) nep->computefunction = fun;
813: if (ctx) nep->functionctx = ctx;
814: if (A) {
815: PetscObjectReference((PetscObject)A);
816: MatDestroy(&nep->function);
817: nep->function = A;
818: }
819: if (B) {
820: PetscObjectReference((PetscObject)B);
821: MatDestroy(&nep->function_pre);
822: nep->function_pre = B;
823: }
824: nep->split = PETSC_FALSE;
825: return(0);
826: }
830: /*@C
831: NEPGetFunction - Returns the Function matrix and optionally the user
832: provided context for evaluating the Function.
834: Not Collective, but Mat object will be parallel if NEP object is
836: Input Parameter:
837: . nep - the nonlinear eigensolver context
839: Output Parameters:
840: + A - location to stash Function matrix (or NULL)
841: . B - location to stash preconditioner matrix (or NULL)
842: . fun - location to put Function function (or NULL)
843: - ctx - location to stash Function context (or NULL)
845: Level: advanced
847: .seealso: NEPSetFunction()
848: @*/
849: PetscErrorCode NEPGetFunction(NEP nep,Mat *A,Mat *B,PetscErrorCode (**fun)(NEP,PetscScalar,Mat*,Mat*,MatStructure*,void*),void **ctx)
850: {
853: if (A) *A = nep->function;
854: if (B) *B = nep->function_pre;
855: if (fun) *fun = nep->computefunction;
856: if (ctx) *ctx = nep->functionctx;
857: return(0);
858: }
862: /*@C
863: NEPSetJacobian - Sets the function to compute Jacobian T'(lambda) as well
864: as the location to store the matrix.
866: Logically Collective on NEP and Mat
868: Input Parameters:
869: + nep - the NEP context
870: . A - Jacobian matrix
871: . jac - Jacobian evaluation routine (if NULL then NEP retains any
872: previously set value)
873: - ctx - [optional] user-defined context for private data for the Jacobian
874: evaluation routine (may be NULL) (if NULL then NEP retains any
875: previously set value)
877: Notes:
878: The routine jac() takes Mat* as the matrix arguments rather than Mat.
879: This allows the Jacobian evaluation routine to replace A with a
880: completely new matrix structure (not just different matrix elements)
881: when appropriate, for instance, if the nonzero structure is changing
882: throughout the global iterations.
884: Level: beginner
886: .seealso: NEPSetFunction(), NEPGetJacobian()
887: @*/
888: PetscErrorCode NEPSetJacobian(NEP nep,Mat A,PetscErrorCode (*jac)(NEP,PetscScalar,Mat*,MatStructure*,void*),void *ctx)
889: {
896: if (jac) nep->computejacobian = jac;
897: if (ctx) nep->jacobianctx = ctx;
898: if (A) {
899: PetscObjectReference((PetscObject)A);
900: MatDestroy(&nep->jacobian);
901: nep->jacobian = A;
902: }
903: nep->split = PETSC_FALSE;
904: return(0);
905: }
909: /*@C
910: NEPGetJacobian - Returns the Jacobian matrix and optionally the user
911: provided context for evaluating the Jacobian.
913: Not Collective, but Mat object will be parallel if NEP object is
915: Input Parameter:
916: . nep - the nonlinear eigensolver context
918: Output Parameters:
919: + A - location to stash Jacobian matrix (or NULL)
920: . jac - location to put Jacobian function (or NULL)
921: - ctx - location to stash Jacobian context (or NULL)
923: Level: advanced
925: .seealso: NEPSetJacobian()
926: @*/
927: PetscErrorCode NEPGetJacobian(NEP nep,Mat *A,PetscErrorCode (**jac)(NEP,PetscScalar,Mat*,MatStructure*,void*),void **ctx)
928: {
931: if (A) *A = nep->jacobian;
932: if (jac) *jac = nep->computejacobian;
933: if (ctx) *ctx = nep->jacobianctx;
934: return(0);
935: }
939: /*@
940: NEPSetSplitOperator - Sets the operator of the nonlinear eigenvalue problem
941: in split form.
943: Collective on NEP, Mat and FN
945: Input Parameters:
946: + nep - the nonlinear eigensolver context
947: . n - number of terms in the split form
948: . A - array of matrices
949: . f - array of functions
950: - str - structure flag for matrices
952: Notes:
953: The nonlinear operator is written as T(lambda) = sum_i A_i*f_i(lambda),
954: for i=1,...,n. The derivative T'(lambda) can be obtained using the
955: derivatives of f_i.
957: The structure flag provides information about A_i's nonzero pattern
958: (see MatStructure enum). If all matrices have the same pattern, then
959: use SAME_NONZERO_PATTERN. If the patterns are different but contained
960: in the pattern of the first one, then use SUBSET_NONZERO_PATTERN.
961: Otherwise use DIFFERENT_NONZERO_PATTERN.
963: This function must be called before NEPSetUp(). If it is called again
964: after NEPSetUp() then the NEP object is reset.
966: Level: intermediate
968: .seealso: NEPGetSplitOperatorTerm(), NEPGetSplitOperatorInfo()
969: @*/
970: PetscErrorCode NEPSetSplitOperator(NEP nep,PetscInt n,Mat A[],FN f[],MatStructure str)
971: {
972: PetscInt i;
978: if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more terms, you have %D",n);
983: if (nep->setupcalled) { NEPReset(nep); }
984: /* clean previously stored information */
985: MatDestroy(&nep->function);
986: MatDestroy(&nep->function_pre);
987: MatDestroy(&nep->jacobian);
988: if (nep->split) {
989: MatDestroyMatrices(nep->nt,&nep->A);
990: for (i=0;i<nep->nt;i++) {
991: FNDestroy(&nep->f[i]);
992: }
993: PetscFree(nep->f);
994: }
995: /* allocate space and copy matrices and functions */
996: PetscMalloc(n*sizeof(Mat),&nep->A);
997: PetscLogObjectMemory(nep,n*sizeof(Mat));
998: for (i=0;i<n;i++) {
1000: PetscObjectReference((PetscObject)A[i]);
1001: nep->A[i] = A[i];
1002: }
1003: PetscMalloc(n*sizeof(FN),&nep->f);
1004: PetscLogObjectMemory(nep,n*sizeof(FN));
1005: for (i=0;i<n;i++) {
1007: PetscObjectReference((PetscObject)f[i]);
1008: nep->f[i] = f[i];
1009: }
1010: nep->nt = n;
1011: nep->mstr = str;
1012: nep->split = PETSC_TRUE;
1013: return(0);
1014: }
1018: /*@
1019: NEPGetSplitOperatorTerm - Gets the matrices and functions associated with
1020: the nonlinear operator in split form.
1022: Not collective, though parallel Mats and FNs are returned if the NEP is parallel
1024: Input Parameter:
1025: + nep - the nonlinear eigensolver context
1026: - k - the index of the requested term (starting in 0)
1028: Output Parameters:
1029: + A - the matrix of the requested term
1030: - f - the function of the requested term
1032: Level: intermediate
1034: .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorInfo()
1035: @*/
1036: PetscErrorCode NEPGetSplitOperatorTerm(NEP nep,PetscInt k,Mat *A,FN *f)
1037: {
1040: if (k<0 || k>=nep->nt) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %d",nep->nt-1);
1041: if (A) *A = nep->A[k];
1042: if (f) *f = nep->f[k];
1043: return(0);
1044: }
1048: /*@
1049: NEPGetSplitOperatorInfo - Returns the number of terms of the split form of
1050: the nonlinear operator, as well as the structure flag for matrices.
1052: Not collective
1054: Input Parameter:
1055: . nep - the nonlinear eigensolver context
1057: Output Parameters:
1058: + n - the number of terms passed in NEPSetSplitOperator()
1059: - str - the matrix structure flag passed in NEPSetSplitOperator()
1061: Level: intermediate
1063: .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorTerm()
1064: @*/
1065: PetscErrorCode NEPGetSplitOperatorInfo(NEP nep,PetscInt *n,MatStructure *str)
1066: {
1069: if (n) *n = nep->nt;
1070: if (str) *str = nep->mstr;
1071: return(0);
1072: }