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) PetscCallBack("PF callback apply to vector", (*pf->ops->applyvec)(pf->data, x, y));
156: else {
157: const PetscScalar *xx;
158: PetscScalar *yy;
160: PetscCall(VecGetLocalSize(x, &n));
161: n = n / pf->dimin;
162: PetscCall(VecGetArrayRead(x, &xx));
163: PetscCall(VecGetArray(y, &yy));
164: PetscCallBack("PF callback apply to array", (*pf->ops->apply)(pf->data, n, xx, yy));
165: PetscCall(VecRestoreArrayRead(x, &xx));
166: PetscCall(VecRestoreArray(y, &yy));
167: }
168: if (nox) PetscCall(VecDestroy(&x));
169: PetscFunctionReturn(PETSC_SUCCESS);
170: }
172: /*@
173: PFApply - Applies the mathematical function to an array of values.
175: Collective
177: Input Parameters:
178: + pf - the function context
179: . n - number of pointwise function evaluations to perform, each pointwise function evaluation
180: is a function of dimin variables and computes dimout variables where dimin and dimout are defined
181: in the call to `PFCreate()`
182: - x - input array
184: Output Parameter:
185: . y - output array
187: Level: beginner
189: .seealso: `PF`, `PFApplyVec()`, `PFCreate()`, `PFDestroy()`, `PFSetType()`, `PFSet()`
190: @*/
191: PetscErrorCode PFApply(PF pf, PetscInt n, const PetscScalar *x, PetscScalar *y)
192: {
193: PetscFunctionBegin;
195: PetscAssertPointer(x, 3);
196: PetscAssertPointer(y, 4);
197: PetscCheck(x != y, PETSC_COMM_SELF, PETSC_ERR_ARG_IDN, "x and y must be different arrays");
199: PetscCallBack("PF callback apply", (*pf->ops->apply)(pf->data, n, x, y));
200: PetscFunctionReturn(PETSC_SUCCESS);
201: }
203: /*@C
204: PFViewFromOptions - View a `PF` based on options set in the options database
206: Collective
208: Input Parameters:
209: + A - the `PF` context
210: . obj - Optional object that provides the prefix used to search the options database
211: - name - command line option
213: Level: intermediate
215: Note:
216: See `PetscObjectViewFromOptions()` for the variety of viewer options available
218: .seealso: `PF`, `PFView`, `PetscObjectViewFromOptions()`, `PFCreate()`
219: @*/
220: PetscErrorCode PFViewFromOptions(PF A, PetscObject obj, const char name[])
221: {
222: PetscFunctionBegin;
224: PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
225: PetscFunctionReturn(PETSC_SUCCESS);
226: }
228: /*@
229: PFView - Prints information about a mathematical function
231: Collective unless `viewer` is `PETSC_VIEWER_STDOUT_SELF`
233: Input Parameters:
234: + pf - the `PF` context
235: - viewer - optional visualization context
237: Level: developer
239: Note:
240: The available visualization contexts include
241: + `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
242: - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
243: output where only the first processor opens
244: the file. All other processors send their
245: data to the first processor to print.
247: The user can open an alternative visualization contexts with
248: `PetscViewerASCIIOpen()` (output to a specified file).
250: .seealso: `PF`, `PetscViewerCreate()`, `PetscViewerASCIIOpen()`
251: @*/
252: PetscErrorCode PFView(PF pf, PetscViewer viewer)
253: {
254: PetscBool iascii;
255: PetscViewerFormat format;
257: PetscFunctionBegin;
259: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pf), &viewer));
261: PetscCheckSameComm(pf, 1, viewer, 2);
263: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
264: if (iascii) {
265: PetscCall(PetscViewerGetFormat(viewer, &format));
266: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pf, viewer));
267: if (pf->ops->view) {
268: PetscCall(PetscViewerASCIIPushTab(viewer));
269: PetscCallBack("PF callback view", (*pf->ops->view)(pf->data, viewer));
270: PetscCall(PetscViewerASCIIPopTab(viewer));
271: }
272: }
273: PetscFunctionReturn(PETSC_SUCCESS);
274: }
276: /*@C
277: PFRegister - Adds a method to the mathematical function package.
279: Not Collective
281: Input Parameters:
282: + sname - name of a new user-defined solver
283: - function - routine to create method context
285: Example Usage:
286: .vb
287: PFRegister("my_function", MyFunctionSetCreate);
288: .ve
290: Then, your solver can be chosen with the procedural interface via
291: $ PFSetType(pf, "my_function")
292: or at runtime via the option
293: $ -pf_type my_function
295: Level: advanced
297: Note:
298: `PFRegister()` may be called multiple times to add several user-defined functions
300: .seealso: `PF`, `PFRegisterAll()`, `PFRegisterDestroy()`
301: @*/
302: PetscErrorCode PFRegister(const char sname[], PetscErrorCode (*function)(PF, void *))
303: {
304: PetscFunctionBegin;
305: PetscCall(PFInitializePackage());
306: PetscCall(PetscFunctionListAdd(&PFList, sname, function));
307: PetscFunctionReturn(PETSC_SUCCESS);
308: }
310: /*@C
311: PFGetType - Gets the `PFType` name (as a string) from the `PF`
312: context.
314: Not Collective
316: Input Parameter:
317: . pf - the function context
319: Output Parameter:
320: . type - name of function
322: Level: intermediate
324: .seealso: `PF`, `PFSetType()`
325: @*/
326: PetscErrorCode PFGetType(PF pf, PFType *type)
327: {
328: PetscFunctionBegin;
330: PetscAssertPointer(type, 2);
331: *type = ((PetscObject)pf)->type_name;
332: PetscFunctionReturn(PETSC_SUCCESS);
333: }
335: /*@C
336: PFSetType - Builds `PF` for a particular function
338: Collective
340: Input Parameters:
341: + pf - the function context.
342: . type - a known method
343: - ctx - optional type dependent context
345: Options Database Key:
346: . -pf_type <type> - Sets PF type
348: Level: intermediate
350: Note:
351: See "petsc/include/petscpf.h" for available methods (for instance, `PFCONSTANT`)
353: .seealso: `PF`, `PFSet()`, `PFRegister()`, `PFCreate()`, `DMDACreatePF()`
354: @*/
355: PetscErrorCode PFSetType(PF pf, PFType type, void *ctx)
356: {
357: PetscBool match;
358: PetscErrorCode (*r)(PF, void *);
360: PetscFunctionBegin;
362: PetscAssertPointer(type, 2);
364: PetscCall(PetscObjectTypeCompare((PetscObject)pf, type, &match));
365: if (match) PetscFunctionReturn(PETSC_SUCCESS);
367: PetscTryTypeMethod(pf, destroy);
368: pf->data = NULL;
370: /* Determine the PFCreateXXX routine for a particular function */
371: PetscCall(PetscFunctionListFind(PFList, type, &r));
372: PetscCheck(r, PetscObjectComm((PetscObject)pf), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested PF type %s", type);
373: pf->ops->destroy = NULL;
374: pf->ops->view = NULL;
375: pf->ops->apply = NULL;
376: pf->ops->applyvec = NULL;
378: /* Call the PFCreateXXX routine for this particular function */
379: PetscCall((*r)(pf, ctx));
381: PetscCall(PetscObjectChangeTypeName((PetscObject)pf, type));
382: PetscFunctionReturn(PETSC_SUCCESS);
383: }
385: /*@
386: PFSetFromOptions - Sets `PF` options from the options database.
388: Collective
390: Input Parameters:
391: . pf - the mathematical function context
393: Level: intermediate
395: Notes:
396: To see all options, run your program with the -help option
397: or consult the users manual.
399: .seealso: `PF`
400: @*/
401: PetscErrorCode PFSetFromOptions(PF pf)
402: {
403: char type[256];
404: PetscBool flg;
406: PetscFunctionBegin;
409: PetscObjectOptionsBegin((PetscObject)pf);
410: PetscCall(PetscOptionsFList("-pf_type", "Type of function", "PFSetType", PFList, NULL, type, 256, &flg));
411: if (flg) PetscCall(PFSetType(pf, type, NULL));
412: PetscTryTypeMethod(pf, setfromoptions, PetscOptionsObject);
414: /* process any options handlers added with PetscObjectAddOptionsHandler() */
415: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)pf, PetscOptionsObject));
416: PetscOptionsEnd();
417: PetscFunctionReturn(PETSC_SUCCESS);
418: }
420: static PetscBool PFPackageInitialized = PETSC_FALSE;
422: /*@C
423: PFFinalizePackage - This function destroys everything in the PETSc `PF` package. It is
424: called from `PetscFinalize()`.
426: Level: developer
428: .seealso: `PF`, `PetscFinalize()`
429: @*/
430: PetscErrorCode PFFinalizePackage(void)
431: {
432: PetscFunctionBegin;
433: PetscCall(PetscFunctionListDestroy(&PFList));
434: PFPackageInitialized = PETSC_FALSE;
435: PFRegisterAllCalled = PETSC_FALSE;
436: PetscFunctionReturn(PETSC_SUCCESS);
437: }
439: /*@C
440: PFInitializePackage - This function initializes everything in the `PF` package. It is called
441: from PetscDLLibraryRegister_petscvec() when using dynamic libraries, and on the first call to `PFCreate()`
442: when using shared or static libraries.
444: Level: developer
446: .seealso: `PF`, `PetscInitialize()`
447: @*/
448: PetscErrorCode PFInitializePackage(void)
449: {
450: char logList[256];
451: PetscBool opt, pkg;
453: PetscFunctionBegin;
454: if (PFPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
455: PFPackageInitialized = PETSC_TRUE;
456: /* Register Classes */
457: PetscCall(PetscClassIdRegister("PointFunction", &PF_CLASSID));
458: /* Register Constructors */
459: PetscCall(PFRegisterAll());
460: /* Process Info */
461: {
462: PetscClassId classids[1];
464: classids[0] = PF_CLASSID;
465: PetscCall(PetscInfoProcessClass("pf", 1, classids));
466: }
467: /* Process summary exclusions */
468: PetscCall(PetscOptionsGetString(NULL, NULL, "-log_exclude", logList, sizeof(logList), &opt));
469: if (opt) {
470: PetscCall(PetscStrInList("pf", logList, ',', &pkg));
471: if (pkg) PetscCall(PetscLogEventExcludeClass(PF_CLASSID));
472: }
473: /* Register package finalizer */
474: PetscCall(PetscRegisterFinalize(PFFinalizePackage));
475: PetscFunctionReturn(PETSC_SUCCESS);
476: }