Actual source code: pf.c
1: /*
2: The PF mathematical functions interface routines, callable by users.
3: */
4: #include <../src/vec/pf/pfimpl.h>
6: PetscClassId PF_CLASSID = 0;
7: PetscFunctionList PFList = NULL; /* list of all registered PD functions */
8: PetscBool PFRegisterAllCalled = PETSC_FALSE;
10: /*@C
11: PFSet - Sets the C/C++/Fortran functions to be used by the PF function
13: Collective
15: Input Parameters:
16: + pf - the function context
17: . apply - function to apply to an array
18: . applyvec - function to apply to a Vec
19: . view - function that prints information about the `PF`
20: . destroy - function to free the private function context
21: - ctx - private function context
23: Level: beginner
25: .seealso: `PF`, `PFCreate()`, `PFDestroy()`, `PFSetType()`, `PFApply()`, `PFApplyVec()`
26: @*/
27: PetscErrorCode PFSet(PF pf, PetscErrorCode (*apply)(void *, PetscInt, const PetscScalar *, PetscScalar *), PetscErrorCode (*applyvec)(void *, Vec, Vec), PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*destroy)(void *), void *ctx)
28: {
29: PetscFunctionBegin;
31: pf->data = ctx;
32: pf->ops->destroy = destroy;
33: pf->ops->apply = apply;
34: pf->ops->applyvec = applyvec;
35: pf->ops->view = view;
36: PetscFunctionReturn(PETSC_SUCCESS);
37: }
39: /*@C
40: PFDestroy - Destroys `PF` context that was created with `PFCreate()`.
42: Collective
44: Input Parameter:
45: . pf - the function context
47: Level: beginner
49: .seealso: `PF`, `PFCreate()`, `PFSet()`, `PFSetType()`
50: @*/
51: PetscErrorCode PFDestroy(PF *pf)
52: {
53: PetscFunctionBegin;
54: if (!*pf) PetscFunctionReturn(PETSC_SUCCESS);
56: if (--((PetscObject)(*pf))->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
58: PetscCall(PFViewFromOptions(*pf, NULL, "-pf_view"));
59: /* if memory was published with SAWs then destroy it */
60: PetscCall(PetscObjectSAWsViewOff((PetscObject)*pf));
62: if ((*pf)->ops->destroy) PetscCall((*(*pf)->ops->destroy)((*pf)->data));
63: PetscCall(PetscHeaderDestroy(pf));
64: PetscFunctionReturn(PETSC_SUCCESS);
65: }
67: /*@C
68: PFCreate - Creates a mathematical function context.
70: Collective
72: Input Parameters:
73: + comm - MPI communicator
74: . dimin - dimension of the space you are mapping from
75: - dimout - dimension of the space you are mapping to
77: Output Parameter:
78: . pf - the function context
80: Level: developer
82: .seealso: `PF`, `PFSet()`, `PFApply()`, `PFDestroy()`, `PFApplyVec()`
83: @*/
84: PetscErrorCode PFCreate(MPI_Comm comm, PetscInt dimin, PetscInt dimout, PF *pf)
85: {
86: PF newpf;
88: PetscFunctionBegin;
89: PetscAssertPointer(pf, 4);
90: *pf = NULL;
91: PetscCall(PFInitializePackage());
93: PetscCall(PetscHeaderCreate(newpf, PF_CLASSID, "PF", "Mathematical functions", "Vec", comm, PFDestroy, PFView));
94: newpf->data = NULL;
95: newpf->ops->destroy = NULL;
96: newpf->ops->apply = NULL;
97: newpf->ops->applyvec = NULL;
98: newpf->ops->view = NULL;
99: newpf->dimin = dimin;
100: newpf->dimout = dimout;
102: *pf = newpf;
103: PetscFunctionReturn(PETSC_SUCCESS);
104: }
106: /* -------------------------------------------------------------------------------*/
108: /*@
109: PFApplyVec - Applies the mathematical function to a vector
111: Collective
113: Input Parameters:
114: + pf - the function context
115: - x - input vector (or `NULL` for the vector (0,1, .... N-1)
117: Output Parameter:
118: . y - output vector
120: Level: beginner
122: .seealso: `PF`, `PFApply()`, `PFCreate()`, `PFDestroy()`, `PFSetType()`, `PFSet()`
123: @*/
124: PetscErrorCode PFApplyVec(PF pf, Vec x, Vec y)
125: {
126: PetscInt i, rstart, rend, n, p;
127: PetscBool nox = PETSC_FALSE;
129: PetscFunctionBegin;
132: if (x) {
134: PetscCheck(x != y, PETSC_COMM_SELF, PETSC_ERR_ARG_IDN, "x and y must be different vectors");
135: } else {
136: PetscScalar *xx;
137: PetscInt lsize;
139: PetscCall(VecGetLocalSize(y, &lsize));
140: lsize = pf->dimin * lsize / pf->dimout;
141: PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)y), lsize, PETSC_DETERMINE, &x));
142: nox = PETSC_TRUE;
143: PetscCall(VecGetOwnershipRange(x, &rstart, &rend));
144: PetscCall(VecGetArray(x, &xx));
145: for (i = rstart; i < rend; i++) xx[i - rstart] = (PetscScalar)i;
146: PetscCall(VecRestoreArray(x, &xx));
147: }
149: PetscCall(VecGetLocalSize(x, &n));
150: PetscCall(VecGetLocalSize(y, &p));
151: PetscCheck((pf->dimin * (n / pf->dimin)) == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local input vector length %" PetscInt_FMT " not divisible by dimin %" PetscInt_FMT " of function", n, pf->dimin);
152: PetscCheck((pf->dimout * (p / pf->dimout)) == p, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local output vector length %" PetscInt_FMT " not divisible by dimout %" PetscInt_FMT " of function", p, pf->dimout);
153: PetscCheck((n / pf->dimin) == (p / pf->dimout), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local vector lengths %" PetscInt_FMT " %" PetscInt_FMT " are wrong for dimin and dimout %" PetscInt_FMT " %" PetscInt_FMT " of function", n, p, pf->dimin, pf->dimout);
155: if (pf->ops->applyvec) PetscCall((*pf->ops->applyvec)(pf->data, x, y));
156: else {
157: PetscScalar *xx, *yy;
159: PetscCall(VecGetLocalSize(x, &n));
160: n = n / pf->dimin;
161: PetscCall(VecGetArray(x, &xx));
162: PetscCall(VecGetArray(y, &yy));
163: PetscCall((*pf->ops->apply)(pf->data, n, xx, yy));
164: PetscCall(VecRestoreArray(x, &xx));
165: PetscCall(VecRestoreArray(y, &yy));
166: }
167: if (nox) PetscCall(VecDestroy(&x));
168: PetscFunctionReturn(PETSC_SUCCESS);
169: }
171: /*@
172: PFApply - Applies the mathematical function to an array of values.
174: Collective
176: Input Parameters:
177: + pf - the function context
178: . n - number of pointwise function evaluations to perform, each pointwise function evaluation
179: is a function of dimin variables and computes dimout variables where dimin and dimout are defined
180: in the call to `PFCreate()`
181: - x - input array
183: Output Parameter:
184: . y - output array
186: Level: beginner
188: .seealso: `PF`, `PFApplyVec()`, `PFCreate()`, `PFDestroy()`, `PFSetType()`, `PFSet()`
189: @*/
190: PetscErrorCode PFApply(PF pf, PetscInt n, const PetscScalar *x, PetscScalar *y)
191: {
192: PetscFunctionBegin;
194: PetscAssertPointer(x, 3);
195: PetscAssertPointer(y, 4);
196: PetscCheck(x != y, PETSC_COMM_SELF, PETSC_ERR_ARG_IDN, "x and y must be different arrays");
198: PetscCall((*pf->ops->apply)(pf->data, n, x, y));
199: PetscFunctionReturn(PETSC_SUCCESS);
200: }
202: /*@C
203: PFViewFromOptions - View a `PF` based on options set in the options database
205: Collective
207: Input Parameters:
208: + A - the `PF` context
209: . obj - Optional object that provides the prefix used to search the options database
210: - name - command line option
212: Level: intermediate
214: Note:
215: See `PetscObjectViewFromOptions()` for the variety of viewer options available
217: .seealso: `PF`, `PFView`, `PetscObjectViewFromOptions()`, `PFCreate()`
218: @*/
219: PetscErrorCode PFViewFromOptions(PF A, PetscObject obj, const char name[])
220: {
221: PetscFunctionBegin;
223: PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
224: PetscFunctionReturn(PETSC_SUCCESS);
225: }
227: /*@
228: PFView - Prints information about a mathematical function
230: Collective unless `viewer` is `PETSC_VIEWER_STDOUT_SELF`
232: Input Parameters:
233: + pf - the `PF` context
234: - viewer - optional visualization context
236: Level: developer
238: Note:
239: The available visualization contexts include
240: + `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
241: - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
242: output where only the first processor opens
243: the file. All other processors send their
244: data to the first processor to print.
246: The user can open an alternative visualization contexts with
247: `PetscViewerASCIIOpen()` (output to a specified file).
249: .seealso: `PF`, `PetscViewerCreate()`, `PetscViewerASCIIOpen()`
250: @*/
251: PetscErrorCode PFView(PF pf, PetscViewer viewer)
252: {
253: PetscBool iascii;
254: PetscViewerFormat format;
256: PetscFunctionBegin;
258: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pf), &viewer));
260: PetscCheckSameComm(pf, 1, viewer, 2);
262: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
263: if (iascii) {
264: PetscCall(PetscViewerGetFormat(viewer, &format));
265: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pf, viewer));
266: if (pf->ops->view) {
267: PetscCall(PetscViewerASCIIPushTab(viewer));
268: PetscCall((*pf->ops->view)(pf->data, viewer));
269: PetscCall(PetscViewerASCIIPopTab(viewer));
270: }
271: }
272: PetscFunctionReturn(PETSC_SUCCESS);
273: }
275: /*@C
276: PFRegister - Adds a method to the mathematical function package.
278: Not Collective
280: Input Parameters:
281: + sname - name of a new user-defined solver
282: - function - routine to create method context
284: Example Usage:
285: .vb
286: PFRegister("my_function", MyFunctionSetCreate);
287: .ve
289: Then, your solver can be chosen with the procedural interface via
290: $ PFSetType(pf, "my_function")
291: or at runtime via the option
292: $ -pf_type my_function
294: Level: advanced
296: Note:
297: `PFRegister()` may be called multiple times to add several user-defined functions
299: .seealso: `PF`, `PFRegisterAll()`, `PFRegisterDestroy()`
300: @*/
301: PetscErrorCode PFRegister(const char sname[], PetscErrorCode (*function)(PF, void *))
302: {
303: PetscFunctionBegin;
304: PetscCall(PFInitializePackage());
305: PetscCall(PetscFunctionListAdd(&PFList, sname, function));
306: PetscFunctionReturn(PETSC_SUCCESS);
307: }
309: /*@C
310: PFGetType - Gets the `PFType` name (as a string) from the `PF`
311: context.
313: Not Collective
315: Input Parameter:
316: . pf - the function context
318: Output Parameter:
319: . type - name of function
321: Level: intermediate
323: .seealso: `PF`, `PFSetType()`
324: @*/
325: PetscErrorCode PFGetType(PF pf, PFType *type)
326: {
327: PetscFunctionBegin;
329: PetscAssertPointer(type, 2);
330: *type = ((PetscObject)pf)->type_name;
331: PetscFunctionReturn(PETSC_SUCCESS);
332: }
334: /*@C
335: PFSetType - Builds `PF` for a particular function
337: Collective
339: Input Parameters:
340: + pf - the function context.
341: . type - a known method
342: - ctx - optional type dependent context
344: Options Database Key:
345: . -pf_type <type> - Sets PF type
347: Level: intermediate
349: Note:
350: See "petsc/include/petscpf.h" for available methods (for instance, `PFCONSTANT`)
352: .seealso: `PF`, `PFSet()`, `PFRegister()`, `PFCreate()`, `DMDACreatePF()`
353: @*/
354: PetscErrorCode PFSetType(PF pf, PFType type, void *ctx)
355: {
356: PetscBool match;
357: PetscErrorCode (*r)(PF, void *);
359: PetscFunctionBegin;
361: PetscAssertPointer(type, 2);
363: PetscCall(PetscObjectTypeCompare((PetscObject)pf, type, &match));
364: if (match) PetscFunctionReturn(PETSC_SUCCESS);
366: PetscTryTypeMethod(pf, destroy);
367: pf->data = NULL;
369: /* Determine the PFCreateXXX routine for a particular function */
370: PetscCall(PetscFunctionListFind(PFList, type, &r));
371: PetscCheck(r, PetscObjectComm((PetscObject)pf), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested PF type %s", type);
372: pf->ops->destroy = NULL;
373: pf->ops->view = NULL;
374: pf->ops->apply = NULL;
375: pf->ops->applyvec = NULL;
377: /* Call the PFCreateXXX routine for this particular function */
378: PetscCall((*r)(pf, ctx));
380: PetscCall(PetscObjectChangeTypeName((PetscObject)pf, type));
381: PetscFunctionReturn(PETSC_SUCCESS);
382: }
384: /*@
385: PFSetFromOptions - Sets `PF` options from the options database.
387: Collective
389: Input Parameters:
390: . pf - the mathematical function context
392: Level: intermediate
394: Notes:
395: To see all options, run your program with the -help option
396: or consult the users manual.
398: .seealso: `PF`
399: @*/
400: PetscErrorCode PFSetFromOptions(PF pf)
401: {
402: char type[256];
403: PetscBool flg;
405: PetscFunctionBegin;
408: PetscObjectOptionsBegin((PetscObject)pf);
409: PetscCall(PetscOptionsFList("-pf_type", "Type of function", "PFSetType", PFList, NULL, type, 256, &flg));
410: if (flg) PetscCall(PFSetType(pf, type, NULL));
411: PetscTryTypeMethod(pf, setfromoptions, PetscOptionsObject);
413: /* process any options handlers added with PetscObjectAddOptionsHandler() */
414: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)pf, PetscOptionsObject));
415: PetscOptionsEnd();
416: PetscFunctionReturn(PETSC_SUCCESS);
417: }
419: static PetscBool PFPackageInitialized = PETSC_FALSE;
421: /*@C
422: PFFinalizePackage - This function destroys everything in the PETSc `PF` package. It is
423: called from `PetscFinalize()`.
425: Level: developer
427: .seealso: `PF`, `PetscFinalize()`
428: @*/
429: PetscErrorCode PFFinalizePackage(void)
430: {
431: PetscFunctionBegin;
432: PetscCall(PetscFunctionListDestroy(&PFList));
433: PFPackageInitialized = PETSC_FALSE;
434: PFRegisterAllCalled = PETSC_FALSE;
435: PetscFunctionReturn(PETSC_SUCCESS);
436: }
438: /*@C
439: PFInitializePackage - This function initializes everything in the `PF` package. It is called
440: from PetscDLLibraryRegister_petscvec() when using dynamic libraries, and on the first call to `PFCreate()`
441: when using shared or static libraries.
443: Level: developer
445: .seealso: `PF`, `PetscInitialize()`
446: @*/
447: PetscErrorCode PFInitializePackage(void)
448: {
449: char logList[256];
450: PetscBool opt, pkg;
452: PetscFunctionBegin;
453: if (PFPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
454: PFPackageInitialized = PETSC_TRUE;
455: /* Register Classes */
456: PetscCall(PetscClassIdRegister("PointFunction", &PF_CLASSID));
457: /* Register Constructors */
458: PetscCall(PFRegisterAll());
459: /* Process Info */
460: {
461: PetscClassId classids[1];
463: classids[0] = PF_CLASSID;
464: PetscCall(PetscInfoProcessClass("pf", 1, classids));
465: }
466: /* Process summary exclusions */
467: PetscCall(PetscOptionsGetString(NULL, NULL, "-log_exclude", logList, sizeof(logList), &opt));
468: if (opt) {
469: PetscCall(PetscStrInList("pf", logList, ',', &pkg));
470: if (pkg) PetscCall(PetscLogEventExcludeClass(PF_CLASSID));
471: }
472: /* Register package finalizer */
473: PetscCall(PetscRegisterFinalize(PFFinalizePackage));
474: PetscFunctionReturn(PETSC_SUCCESS);
475: }