Actual source code: dmts.c
1: #include <petsc/private/tsimpl.h>
2: #include <petsc/private/dmimpl.h>
4: static PetscErrorCode DMTSUnsetRHSFunctionContext_DMTS(DMTS tsdm)
5: {
6: PetscFunctionBegin;
7: PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", NULL));
8: tsdm->rhsfunctionctxcontainer = NULL;
9: PetscFunctionReturn(PETSC_SUCCESS);
10: }
12: static PetscErrorCode DMTSUnsetRHSJacobianContext_DMTS(DMTS tsdm)
13: {
14: PetscFunctionBegin;
15: PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", NULL));
16: tsdm->rhsjacobianctxcontainer = NULL;
17: PetscFunctionReturn(PETSC_SUCCESS);
18: }
20: static PetscErrorCode DMTSUnsetIFunctionContext_DMTS(DMTS tsdm)
21: {
22: PetscFunctionBegin;
23: PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", NULL));
24: tsdm->ifunctionctxcontainer = NULL;
25: PetscFunctionReturn(PETSC_SUCCESS);
26: }
28: static PetscErrorCode DMTSUnsetIJacobianContext_DMTS(DMTS tsdm)
29: {
30: PetscFunctionBegin;
31: PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", NULL));
32: tsdm->ijacobianctxcontainer = NULL;
33: PetscFunctionReturn(PETSC_SUCCESS);
34: }
36: static PetscErrorCode DMTSUnsetI2FunctionContext_DMTS(DMTS tsdm)
37: {
38: PetscFunctionBegin;
39: PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", NULL));
40: tsdm->i2functionctxcontainer = NULL;
41: PetscFunctionReturn(PETSC_SUCCESS);
42: }
44: static PetscErrorCode DMTSUnsetI2JacobianContext_DMTS(DMTS tsdm)
45: {
46: PetscFunctionBegin;
47: PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", NULL));
48: tsdm->i2jacobianctxcontainer = NULL;
49: PetscFunctionReturn(PETSC_SUCCESS);
50: }
52: static PetscErrorCode DMTSDestroy(DMTS *kdm)
53: {
54: PetscFunctionBegin;
55: if (!*kdm) PetscFunctionReturn(PETSC_SUCCESS);
57: if (--((PetscObject)*kdm)->refct > 0) {
58: *kdm = NULL;
59: PetscFunctionReturn(PETSC_SUCCESS);
60: }
61: PetscCall(DMTSUnsetRHSFunctionContext_DMTS(*kdm));
62: PetscCall(DMTSUnsetRHSJacobianContext_DMTS(*kdm));
63: PetscCall(DMTSUnsetIFunctionContext_DMTS(*kdm));
64: PetscCall(DMTSUnsetIJacobianContext_DMTS(*kdm));
65: PetscCall(DMTSUnsetI2FunctionContext_DMTS(*kdm));
66: PetscCall(DMTSUnsetI2JacobianContext_DMTS(*kdm));
67: PetscTryTypeMethod(*kdm, destroy);
68: PetscCall(PetscHeaderDestroy(kdm));
69: PetscFunctionReturn(PETSC_SUCCESS);
70: }
72: PetscErrorCode DMTSLoad(DMTS kdm, PetscViewer viewer)
73: {
74: PetscFunctionBegin;
75: PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunction, 1, NULL, PETSC_FUNCTION));
76: PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunctionview, 1, NULL, PETSC_FUNCTION));
77: PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunctionload, 1, NULL, PETSC_FUNCTION));
78: if (kdm->ops->ifunctionload) {
79: void *ctx;
81: PetscCall(PetscContainerGetPointer(kdm->ifunctionctxcontainer, &ctx));
82: PetscCall((*kdm->ops->ifunctionload)(&ctx, viewer));
83: }
84: PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobian, 1, NULL, PETSC_FUNCTION));
85: PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobianview, 1, NULL, PETSC_FUNCTION));
86: PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobianload, 1, NULL, PETSC_FUNCTION));
87: if (kdm->ops->ijacobianload) {
88: void *ctx;
90: PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer, &ctx));
91: PetscCall((*kdm->ops->ijacobianload)(&ctx, viewer));
92: }
93: PetscFunctionReturn(PETSC_SUCCESS);
94: }
96: PetscErrorCode DMTSView(DMTS kdm, PetscViewer viewer)
97: {
98: PetscBool isascii, isbinary;
100: PetscFunctionBegin;
101: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
102: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
103: if (isascii) {
104: #if defined(PETSC_SERIALIZE_FUNCTIONS)
105: const char *fname;
107: PetscCall(PetscFPTFind(kdm->ops->ifunction, &fname));
108: if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, " IFunction used by TS: %s\n", fname));
109: PetscCall(PetscFPTFind(kdm->ops->ijacobian, &fname));
110: if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, " IJacobian function used by TS: %s\n", fname));
111: #endif
112: } else if (isbinary) {
113: struct {
114: TSIFunctionFn *ifunction;
115: } funcstruct;
116: struct {
117: PetscErrorCode (*ifunctionview)(void *, PetscViewer);
118: } funcviewstruct;
119: struct {
120: PetscErrorCode (*ifunctionload)(void **, PetscViewer);
121: } funcloadstruct;
122: struct {
123: TSIJacobianFn *ijacobian;
124: } jacstruct;
125: struct {
126: PetscErrorCode (*ijacobianview)(void *, PetscViewer);
127: } jacviewstruct;
128: struct {
129: PetscErrorCode (*ijacobianload)(void **, PetscViewer);
130: } jacloadstruct;
132: funcstruct.ifunction = kdm->ops->ifunction;
133: funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
134: funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
135: PetscCall(PetscViewerBinaryWrite(viewer, &funcstruct, 1, PETSC_FUNCTION));
136: PetscCall(PetscViewerBinaryWrite(viewer, &funcviewstruct, 1, PETSC_FUNCTION));
137: PetscCall(PetscViewerBinaryWrite(viewer, &funcloadstruct, 1, PETSC_FUNCTION));
138: if (kdm->ops->ifunctionview) {
139: void *ctx;
141: PetscCall(PetscContainerGetPointer(kdm->ifunctionctxcontainer, &ctx));
142: PetscCall((*kdm->ops->ifunctionview)(ctx, viewer));
143: }
144: jacstruct.ijacobian = kdm->ops->ijacobian;
145: jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
146: jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
147: PetscCall(PetscViewerBinaryWrite(viewer, &jacstruct, 1, PETSC_FUNCTION));
148: PetscCall(PetscViewerBinaryWrite(viewer, &jacviewstruct, 1, PETSC_FUNCTION));
149: PetscCall(PetscViewerBinaryWrite(viewer, &jacloadstruct, 1, PETSC_FUNCTION));
150: if (kdm->ops->ijacobianview) {
151: void *ctx;
153: PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer, &ctx));
154: PetscCall((*kdm->ops->ijacobianview)(ctx, viewer));
155: }
156: }
157: PetscFunctionReturn(PETSC_SUCCESS);
158: }
160: static PetscErrorCode DMTSCreate(MPI_Comm comm, DMTS *kdm)
161: {
162: PetscFunctionBegin;
163: PetscCall(TSInitializePackage());
164: PetscCall(PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView));
165: PetscFunctionReturn(PETSC_SUCCESS);
166: }
168: /* Attaches the DMTS to the coarse level.
169: * Under what conditions should we copy versus duplicate?
170: */
171: static PetscErrorCode DMCoarsenHook_DMTS(DM dm, DM dmc, void *ctx)
172: {
173: PetscFunctionBegin;
174: PetscCall(DMCopyDMTS(dm, dmc));
175: PetscFunctionReturn(PETSC_SUCCESS);
176: }
178: /* This could restrict auxiliary information to the coarse level.
179: */
180: static PetscErrorCode DMRestrictHook_DMTS(DM dm, Mat Restrict, Vec rscale, Mat Inject, DM dmc, void *ctx)
181: {
182: PetscFunctionBegin;
183: PetscFunctionReturn(PETSC_SUCCESS);
184: }
186: static PetscErrorCode DMSubDomainHook_DMTS(DM dm, DM subdm, void *ctx)
187: {
188: PetscFunctionBegin;
189: PetscCall(DMCopyDMTS(dm, subdm));
190: PetscFunctionReturn(PETSC_SUCCESS);
191: }
193: /* This could restrict auxiliary information to the coarse level.
194: */
195: static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
196: {
197: PetscFunctionBegin;
198: PetscFunctionReturn(PETSC_SUCCESS);
199: }
201: /*@C
202: DMTSCopy - copies the information in a `DMTS` to another `DMTS`
204: Not Collective
206: Input Parameters:
207: + kdm - Original `DMTS`
208: - nkdm - `DMTS` to receive the data, should have been created with `DMTSCreate()`
210: Level: developer
212: .seealso: [](ch_ts), `DMTSCreate()`, `DMTSDestroy()`
213: @*/
214: PetscErrorCode DMTSCopy(DMTS kdm, DMTS nkdm)
215: {
216: PetscFunctionBegin;
219: nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
220: nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
221: nkdm->ops->ifunction = kdm->ops->ifunction;
222: nkdm->ops->ijacobian = kdm->ops->ijacobian;
223: nkdm->ops->i2function = kdm->ops->i2function;
224: nkdm->ops->i2jacobian = kdm->ops->i2jacobian;
225: nkdm->ops->solution = kdm->ops->solution;
226: nkdm->ops->destroy = kdm->ops->destroy;
227: nkdm->ops->duplicate = kdm->ops->duplicate;
229: nkdm->solutionctx = kdm->solutionctx;
230: nkdm->rhsfunctionctxcontainer = kdm->rhsfunctionctxcontainer;
231: nkdm->rhsjacobianctxcontainer = kdm->rhsjacobianctxcontainer;
232: nkdm->ifunctionctxcontainer = kdm->ifunctionctxcontainer;
233: nkdm->ijacobianctxcontainer = kdm->ijacobianctxcontainer;
234: nkdm->i2functionctxcontainer = kdm->i2functionctxcontainer;
235: nkdm->i2jacobianctxcontainer = kdm->i2jacobianctxcontainer;
236: if (nkdm->rhsfunctionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "rhs function ctx", (PetscObject)nkdm->rhsfunctionctxcontainer));
237: if (nkdm->rhsjacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "rhs jacobian ctx", (PetscObject)nkdm->rhsjacobianctxcontainer));
238: if (nkdm->ifunctionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "ifunction ctx", (PetscObject)nkdm->ifunctionctxcontainer));
239: if (nkdm->ijacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "ijacobian ctx", (PetscObject)nkdm->ijacobianctxcontainer));
240: if (nkdm->i2functionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "i2function ctx", (PetscObject)nkdm->i2functionctxcontainer));
241: if (nkdm->i2jacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "i2jacobian ctx", (PetscObject)nkdm->i2jacobianctxcontainer));
243: nkdm->data = kdm->data;
245: /*
246: nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
247: nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
248: nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
249: */
251: /* implementation specific copy hooks */
252: PetscTryTypeMethod(kdm, duplicate, nkdm);
253: PetscFunctionReturn(PETSC_SUCCESS);
254: }
256: /*@C
257: DMGetDMTS - get read-only private `DMTS` context from a `DM`
259: Not Collective
261: Input Parameter:
262: . dm - `DM` to be used with `TS`
264: Output Parameter:
265: . tsdm - private `DMTS` context
267: Level: developer
269: Notes:
270: Use `DMGetDMTSWrite()` if write access is needed. The `DMTSSetXXX()` API should be used wherever possible.
272: .seealso: [](ch_ts), `DMTS`, `DMGetDMTSWrite()`
273: @*/
274: PetscErrorCode DMGetDMTS(DM dm, DMTS *tsdm)
275: {
276: PetscFunctionBegin;
278: *tsdm = (DMTS)dm->dmts;
279: if (!*tsdm) {
280: PetscCall(PetscInfo(dm, "Creating new DMTS\n"));
281: PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm), tsdm));
282: dm->dmts = (PetscObject)*tsdm;
283: (*tsdm)->originaldm = dm;
284: PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_DMTS, DMRestrictHook_DMTS, NULL));
285: PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_DMTS, DMSubDomainRestrictHook_DMTS, NULL));
286: }
287: PetscFunctionReturn(PETSC_SUCCESS);
288: }
290: /*@C
291: DMGetDMTSWrite - get write access to private `DMTS` context from a `DM`
293: Not Collective
295: Input Parameter:
296: . dm - `DM` to be used with `TS`
298: Output Parameter:
299: . tsdm - private `DMTS` context
301: Level: developer
303: .seealso: [](ch_ts), `DMTS`, `DMGetDMTS()`
304: @*/
305: PetscErrorCode DMGetDMTSWrite(DM dm, DMTS *tsdm)
306: {
307: DMTS sdm;
309: PetscFunctionBegin;
311: PetscCall(DMGetDMTS(dm, &sdm));
312: PetscCheck(sdm->originaldm, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DMTS has a NULL originaldm");
313: if (sdm->originaldm != dm) { /* Copy on write */
314: DMTS oldsdm = sdm;
315: PetscCall(PetscInfo(dm, "Copying DMTS due to write\n"));
316: PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm), &sdm));
317: PetscCall(DMTSCopy(oldsdm, sdm));
318: PetscCall(DMTSDestroy((DMTS *)&dm->dmts));
319: dm->dmts = (PetscObject)sdm;
320: sdm->originaldm = dm;
321: }
322: *tsdm = sdm;
323: PetscFunctionReturn(PETSC_SUCCESS);
324: }
326: /*@C
327: DMCopyDMTS - copies a `DMTS` context to a new `DM`
329: Logically Collective
331: Input Parameters:
332: + dmsrc - `DM` to obtain context from
333: - dmdest - `DM` to add context to
335: Level: developer
337: Note:
338: The context is copied by reference. This function does not ensure that a context exists.
340: .seealso: [](ch_ts), `DMTS`, `DMGetDMTS()`, `TSSetDM()`
341: @*/
342: PetscErrorCode DMCopyDMTS(DM dmsrc, DM dmdest)
343: {
344: PetscFunctionBegin;
347: PetscCall(DMTSDestroy((DMTS *)&dmdest->dmts));
348: dmdest->dmts = dmsrc->dmts;
349: PetscCall(PetscObjectReference(dmdest->dmts));
350: PetscCall(DMCoarsenHookAdd(dmdest, DMCoarsenHook_DMTS, DMRestrictHook_DMTS, NULL));
351: PetscCall(DMSubDomainHookAdd(dmdest, DMSubDomainHook_DMTS, DMSubDomainRestrictHook_DMTS, NULL));
352: PetscFunctionReturn(PETSC_SUCCESS);
353: }
355: /*@C
356: DMTSSetIFunction - set `TS` implicit function evaluation function into a `DMTS`
358: Not Collective
360: Input Parameters:
361: + dm - `DM` to be used with `TS`
362: . func - function evaluating f(t,u,u_t)
363: - ctx - context for residual evaluation
365: Level: developer
367: Note:
368: `TSSetIFunction()` is normally used, but it calls this function internally because the user context is actually
369: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
370: not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
372: .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `TSIFunctionFn`
373: @*/
374: PetscErrorCode DMTSSetIFunction(DM dm, TSIFunctionFn *func, void *ctx)
375: {
376: DMTS tsdm;
378: PetscFunctionBegin;
380: PetscCall(DMGetDMTSWrite(dm, &tsdm));
381: if (func) tsdm->ops->ifunction = func;
382: if (ctx) {
383: PetscContainer ctxcontainer;
384: PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
385: PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
386: PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", (PetscObject)ctxcontainer));
387: tsdm->ifunctionctxcontainer = ctxcontainer;
388: PetscCall(PetscContainerDestroy(&ctxcontainer));
389: }
390: PetscFunctionReturn(PETSC_SUCCESS);
391: }
393: /*@C
394: DMTSSetIFunctionContextDestroy - set `TS` implicit evaluation context destroy function into a `DMTS`
396: Not Collective
398: Input Parameters:
399: + dm - `DM` to be used with `TS`
400: - f - implicit evaluation context destroy function
402: Level: developer
404: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetIFunction()`, `TSSetIFunction()`
405: @*/
406: PetscErrorCode DMTSSetIFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
407: {
408: DMTS tsdm;
410: PetscFunctionBegin;
412: PetscCall(DMGetDMTSWrite(dm, &tsdm));
413: if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ifunctionctxcontainer, f));
414: PetscFunctionReturn(PETSC_SUCCESS);
415: }
417: PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM dm)
418: {
419: DMTS tsdm;
421: PetscFunctionBegin;
423: PetscCall(DMGetDMTSWrite(dm, &tsdm));
424: PetscCall(DMTSUnsetIFunctionContext_DMTS(tsdm));
425: PetscFunctionReturn(PETSC_SUCCESS);
426: }
428: /*@C
429: DMTSGetIFunction - get `TS` implicit residual evaluation function from a `DMTS`
431: Not Collective
433: Input Parameter:
434: . dm - `DM` to be used with `TS`
436: Output Parameters:
437: + func - function evaluation function, for calling sequence see `TSIFunctionFn`
438: - ctx - context for residual evaluation
440: Level: developer
442: Note:
443: `TSGetIFunction()` is normally used, but it calls this function internally because the user context is actually
444: associated with the `DM`.
446: .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `DMTSSetIFunction()`, `TSIFunctionFn`
447: @*/
448: PetscErrorCode DMTSGetIFunction(DM dm, TSIFunctionFn **func, void **ctx)
449: {
450: DMTS tsdm;
452: PetscFunctionBegin;
454: PetscCall(DMGetDMTS(dm, &tsdm));
455: if (func) *func = tsdm->ops->ifunction;
456: if (ctx) {
457: if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ifunctionctxcontainer, ctx));
458: else *ctx = NULL;
459: }
460: PetscFunctionReturn(PETSC_SUCCESS);
461: }
463: /*@C
464: DMTSSetI2Function - set `TS` implicit function evaluation function for 2nd order systems into a `TSDM`
466: Not Collective
468: Input Parameters:
469: + dm - `DM` to be used with `TS`
470: . fun - function evaluation routine
471: - ctx - context for residual evaluation
473: Level: developer
475: Note:
476: `TSSetI2Function()` is normally used, but it calls this function internally because the user context is actually
477: associated with the `DM`.
479: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSSetI2Function()`
480: @*/
481: PetscErrorCode DMTSSetI2Function(DM dm, TSI2FunctionFn *fun, void *ctx)
482: {
483: DMTS tsdm;
485: PetscFunctionBegin;
487: PetscCall(DMGetDMTSWrite(dm, &tsdm));
488: if (fun) tsdm->ops->i2function = fun;
489: if (ctx) {
490: PetscContainer ctxcontainer;
491: PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
492: PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
493: PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", (PetscObject)ctxcontainer));
494: tsdm->i2functionctxcontainer = ctxcontainer;
495: PetscCall(PetscContainerDestroy(&ctxcontainer));
496: }
497: PetscFunctionReturn(PETSC_SUCCESS);
498: }
500: /*@C
501: DMTSSetI2FunctionContextDestroy - set `TS` implicit evaluation for 2nd order systems context destroy into a `DMTS`
503: Not Collective
505: Input Parameters:
506: + dm - `DM` to be used with `TS`
507: - f - implicit evaluation context destroy function
509: Level: developer
511: Note:
512: `TSSetI2FunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
513: associated with the `DM`.
515: .seealso: [](ch_ts), `DMTS`, `TSSetI2FunctionContextDestroy()`, `DMTSSetI2Function()`, `TSSetI2Function()`
516: @*/
517: PetscErrorCode DMTSSetI2FunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
518: {
519: DMTS tsdm;
521: PetscFunctionBegin;
523: PetscCall(DMGetDMTSWrite(dm, &tsdm));
524: if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2functionctxcontainer, f));
525: PetscFunctionReturn(PETSC_SUCCESS);
526: }
528: PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM dm)
529: {
530: DMTS tsdm;
532: PetscFunctionBegin;
534: PetscCall(DMGetDMTSWrite(dm, &tsdm));
535: PetscCall(DMTSUnsetI2FunctionContext_DMTS(tsdm));
536: PetscFunctionReturn(PETSC_SUCCESS);
537: }
539: /*@C
540: DMTSGetI2Function - get `TS` implicit residual evaluation function for 2nd order systems from a `DMTS`
542: Not Collective
544: Input Parameter:
545: . dm - `DM` to be used with `TS`
547: Output Parameters:
548: + fun - function evaluation function, for calling sequence see `TSSetI2Function()`
549: - ctx - context for residual evaluation
551: Level: developer
553: Note:
554: `TSGetI2Function()` is normally used, but it calls this function internally because the user context is actually
555: associated with the `DM`.
557: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetI2Function()`, `TSGetI2Function()`
558: @*/
559: PetscErrorCode DMTSGetI2Function(DM dm, TSI2FunctionFn **fun, void **ctx)
560: {
561: DMTS tsdm;
563: PetscFunctionBegin;
565: PetscCall(DMGetDMTS(dm, &tsdm));
566: if (fun) *fun = tsdm->ops->i2function;
567: if (ctx) {
568: if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2functionctxcontainer, ctx));
569: else *ctx = NULL;
570: }
571: PetscFunctionReturn(PETSC_SUCCESS);
572: }
574: /*@C
575: DMTSSetI2Jacobian - set `TS` implicit Jacobian evaluation function for 2nd order systems from a `DMTS`
577: Not Collective
579: Input Parameters:
580: + dm - `DM` to be used with `TS`
581: . jac - Jacobian evaluation routine
582: - ctx - context for Jacobian evaluation
584: Level: developer
586: Note:
587: `TSSetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
588: associated with the `DM`.
590: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSI2JacobianFn`, `TSSetI2Jacobian()`
591: @*/
592: PetscErrorCode DMTSSetI2Jacobian(DM dm, TSI2JacobianFn *jac, void *ctx)
593: {
594: DMTS tsdm;
596: PetscFunctionBegin;
598: PetscCall(DMGetDMTSWrite(dm, &tsdm));
599: if (jac) tsdm->ops->i2jacobian = jac;
600: if (ctx) {
601: PetscContainer ctxcontainer;
602: PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
603: PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
604: PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", (PetscObject)ctxcontainer));
605: tsdm->i2jacobianctxcontainer = ctxcontainer;
606: PetscCall(PetscContainerDestroy(&ctxcontainer));
607: }
608: PetscFunctionReturn(PETSC_SUCCESS);
609: }
611: /*@C
612: DMTSSetI2JacobianContextDestroy - set `TS` implicit Jacobian evaluation for 2nd order systems context destroy function into a `DMTS`
614: Not Collective
616: Input Parameters:
617: + dm - `DM` to be used with `TS`
618: - f - implicit Jacobian evaluation context destroy function
620: Level: developer
622: Note:
623: Normally `TSSetI2JacobianContextDestroy()` is used
625: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSSetI2JacobianContextDestroy()`, `DMTSSetI2Jacobian()`, `TSSetI2Jacobian()`
626: @*/
627: PetscErrorCode DMTSSetI2JacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
628: {
629: DMTS tsdm;
631: PetscFunctionBegin;
633: PetscCall(DMGetDMTSWrite(dm, &tsdm));
634: if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2jacobianctxcontainer, f));
635: PetscFunctionReturn(PETSC_SUCCESS);
636: }
638: PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM dm)
639: {
640: DMTS tsdm;
642: PetscFunctionBegin;
644: PetscCall(DMGetDMTSWrite(dm, &tsdm));
645: PetscCall(DMTSUnsetI2JacobianContext_DMTS(tsdm));
646: PetscFunctionReturn(PETSC_SUCCESS);
647: }
649: /*@C
650: DMTSGetI2Jacobian - get `TS` implicit Jacobian evaluation function for 2nd order systems from a `DMTS`
652: Not Collective
654: Input Parameter:
655: . dm - `DM` to be used with `TS`
657: Output Parameters:
658: + jac - Jacobian evaluation function, for calling sequence see `TSI2JacobianFn`
659: - ctx - context for Jacobian evaluation
661: Level: developer
663: Note:
664: `TSGetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
665: associated with the `DM`.
667: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()`, `TSI2JacobianFn`
668: @*/
669: PetscErrorCode DMTSGetI2Jacobian(DM dm, TSI2JacobianFn **jac, void **ctx)
670: {
671: DMTS tsdm;
673: PetscFunctionBegin;
675: PetscCall(DMGetDMTS(dm, &tsdm));
676: if (jac) *jac = tsdm->ops->i2jacobian;
677: if (ctx) {
678: if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2jacobianctxcontainer, ctx));
679: else *ctx = NULL;
680: }
681: PetscFunctionReturn(PETSC_SUCCESS);
682: }
684: /*@C
685: DMTSSetRHSFunction - set `TS` explicit residual evaluation function into a `DMTS`
687: Not Collective
689: Input Parameters:
690: + dm - `DM` to be used with `TS`
691: . func - RHS function evaluation routine, see `TSRHSFunctionFn` for the calling sequence
692: - ctx - context for residual evaluation
694: Level: developer
696: Note:
697: `TSSetRHSFunction()` is normally used, but it calls this function internally because the user context is actually
698: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
699: not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
701: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSRHSFunctionFn`
702: @*/
703: PetscErrorCode DMTSSetRHSFunction(DM dm, TSRHSFunctionFn *func, void *ctx)
704: {
705: DMTS tsdm;
707: PetscFunctionBegin;
709: PetscCall(DMGetDMTSWrite(dm, &tsdm));
710: if (func) tsdm->ops->rhsfunction = func;
711: if (ctx) {
712: PetscContainer ctxcontainer;
713: PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
714: PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
715: PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", (PetscObject)ctxcontainer));
716: tsdm->rhsfunctionctxcontainer = ctxcontainer;
717: PetscCall(PetscContainerDestroy(&ctxcontainer));
718: }
719: PetscFunctionReturn(PETSC_SUCCESS);
720: }
722: /*@C
723: DMTSSetRHSFunctionContextDestroy - set `TS` explicit residual evaluation context destroy function into a `DMTS`
725: Not Collective
727: Input Parameters:
728: + dm - `DM` to be used with `TS`
729: - f - explicit evaluation context destroy function
731: Level: developer
733: Note:
734: `TSSetRHSFunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
735: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
736: not.
738: Developer Notes:
739: If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
741: .seealso: [](ch_ts), `DMTS`, `TSSetRHSFunctionContextDestroy()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
742: @*/
743: PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
744: {
745: DMTS tsdm;
747: PetscFunctionBegin;
749: PetscCall(DMGetDMTSWrite(dm, &tsdm));
750: if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsfunctionctxcontainer, f));
751: PetscFunctionReturn(PETSC_SUCCESS);
752: }
754: PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM dm)
755: {
756: DMTS tsdm;
758: PetscFunctionBegin;
760: PetscCall(DMGetDMTSWrite(dm, &tsdm));
761: PetscCall(DMTSUnsetRHSFunctionContext_DMTS(tsdm));
762: tsdm->rhsfunctionctxcontainer = NULL;
763: PetscFunctionReturn(PETSC_SUCCESS);
764: }
766: /*@C
767: DMTSSetTransientVariable - sets function to transform from state to transient variables into a `DMTS`
769: Logically Collective
771: Input Parameters:
772: + dm - `DM` to be used with `TS`
773: . tvar - a function that transforms to transient variables, see `TSTransientVariableFn` for the calling sequence
774: - ctx - a context for tvar
776: Level: developer
778: Notes:
779: Normally `TSSetTransientVariable()` is used
781: This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., `TSBDF`)
782: can be conservative. In this context, primitive variables P are used to model the state (e.g., because they lead to
783: well-conditioned formulations even in limiting cases such as low-Mach or zero porosity). The transient variable is
784: C(P), specified by calling this function. An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
785: evaluated via the chain rule, as in
787: $$
788: dF/dP + shift * dF/dCdot dC/dP.
789: $$
791: .seealso: [](ch_ts), `DMTS`, `TS`, `TSBDF`, `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()`, `TSTransientVariableFn`
792: @*/
793: PetscErrorCode DMTSSetTransientVariable(DM dm, TSTransientVariableFn *tvar, void *ctx)
794: {
795: DMTS dmts;
797: PetscFunctionBegin;
799: PetscCall(DMGetDMTSWrite(dm, &dmts));
800: dmts->ops->transientvar = tvar;
801: dmts->transientvarctx = ctx;
802: PetscFunctionReturn(PETSC_SUCCESS);
803: }
805: /*@C
806: DMTSGetTransientVariable - gets function to transform from state to transient variables set with `DMTSSetTransientVariable()` from a `TSDM`
808: Logically Collective
810: Input Parameter:
811: . dm - `DM` to be used with `TS`
813: Output Parameters:
814: + tvar - a function that transforms to transient variables, see `TSTransientVariableFn` for the calling sequence
815: - ctx - a context for tvar
817: Level: developer
819: Note:
820: Normally `TSSetTransientVariable()` is used
822: .seealso: [](ch_ts), `DMTS`, `DM`, `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()`, `TSTransientVariableFn`
823: @*/
824: PetscErrorCode DMTSGetTransientVariable(DM dm, TSTransientVariableFn **tvar, void *ctx)
825: {
826: DMTS dmts;
828: PetscFunctionBegin;
830: PetscCall(DMGetDMTS(dm, &dmts));
831: if (tvar) *tvar = dmts->ops->transientvar;
832: if (ctx) *(void **)ctx = dmts->transientvarctx;
833: PetscFunctionReturn(PETSC_SUCCESS);
834: }
836: /*@C
837: DMTSGetSolutionFunction - gets the `TS` solution evaluation function from a `DMTS`
839: Not Collective
841: Input Parameter:
842: . dm - `DM` to be used with `TS`
844: Output Parameters:
845: + func - solution function evaluation function, for calling sequence see `TSSolutionFn`
846: - ctx - context for solution evaluation
848: Level: developer
850: .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `DMTSSetSolutionFunction()`, `TSSolutionFn`
851: @*/
852: PetscErrorCode DMTSGetSolutionFunction(DM dm, TSSolutionFn **func, void **ctx)
853: {
854: DMTS tsdm;
856: PetscFunctionBegin;
858: PetscCall(DMGetDMTS(dm, &tsdm));
859: if (func) *func = tsdm->ops->solution;
860: if (ctx) *ctx = tsdm->solutionctx;
861: PetscFunctionReturn(PETSC_SUCCESS);
862: }
864: /*@C
865: DMTSSetSolutionFunction - set `TS` solution evaluation function into a `DMTS`
867: Not Collective
869: Input Parameters:
870: + dm - `DM` to be used with `TS`
871: . func - solution function evaluation routine, for calling sequence see `TSSolutionFn`
872: - ctx - context for solution evaluation
874: Level: developer
876: Note:
877: `TSSetSolutionFunction()` is normally used, but it calls this function internally because the user context is actually
878: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
879: not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
881: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSGetSolutionFunction()`, `TSSolutionFn`
882: @*/
883: PetscErrorCode DMTSSetSolutionFunction(DM dm, TSSolutionFn *func, void *ctx)
884: {
885: DMTS tsdm;
887: PetscFunctionBegin;
889: PetscCall(DMGetDMTSWrite(dm, &tsdm));
890: if (func) tsdm->ops->solution = func;
891: if (ctx) tsdm->solutionctx = ctx;
892: PetscFunctionReturn(PETSC_SUCCESS);
893: }
895: /*@C
896: DMTSSetForcingFunction - set `TS` forcing function evaluation function into a `DMTS`
898: Not Collective
900: Input Parameters:
901: + dm - `DM` to be used with `TS`
902: . func - forcing function evaluation routine, for calling sequence see `TSForcingFn`
903: - ctx - context for solution evaluation
905: Level: developer
907: Note:
908: `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually
909: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
910: not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
912: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSForcingFn`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
913: @*/
914: PetscErrorCode DMTSSetForcingFunction(DM dm, TSForcingFn *func, void *ctx)
915: {
916: DMTS tsdm;
918: PetscFunctionBegin;
920: PetscCall(DMGetDMTSWrite(dm, &tsdm));
921: if (func) tsdm->ops->forcing = func;
922: if (ctx) tsdm->forcingctx = ctx;
923: PetscFunctionReturn(PETSC_SUCCESS);
924: }
926: /*@C
927: DMTSGetForcingFunction - get `TS` forcing function evaluation function from a `DMTS`
929: Not Collective
931: Input Parameter:
932: . dm - `DM` to be used with `TS`
934: Output Parameters:
935: + f - forcing function evaluation function; see `TSForcingFn` for the calling sequence
936: - ctx - context for solution evaluation
938: Level: developer
940: Note:
941: `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually
942: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
943: not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.
945: .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `TSSetForcingFunction()`, `TSForcingFn`
946: @*/
947: PetscErrorCode DMTSGetForcingFunction(DM dm, TSForcingFn **f, void **ctx)
948: {
949: DMTS tsdm;
951: PetscFunctionBegin;
953: PetscCall(DMGetDMTSWrite(dm, &tsdm));
954: if (f) *f = tsdm->ops->forcing;
955: if (ctx) *ctx = tsdm->forcingctx;
956: PetscFunctionReturn(PETSC_SUCCESS);
957: }
959: /*@C
960: DMTSGetRHSFunction - get `TS` explicit residual evaluation function from a `DMTS`
962: Not Collective
964: Input Parameter:
965: . dm - `DM` to be used with `TS`
967: Output Parameters:
968: + func - residual evaluation function, for calling sequence see `TSRHSFunctionFn`
969: - ctx - context for residual evaluation
971: Level: developer
973: Note:
974: `TSGetRHSFunction()` is normally used, but it calls this function internally because the user context is actually
975: associated with the DM.
977: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSRHSFunctionFn`, `TSGetRHSFunction()`
978: @*/
979: PetscErrorCode DMTSGetRHSFunction(DM dm, TSRHSFunctionFn **func, void **ctx)
980: {
981: DMTS tsdm;
983: PetscFunctionBegin;
985: PetscCall(DMGetDMTS(dm, &tsdm));
986: if (func) *func = tsdm->ops->rhsfunction;
987: if (ctx) {
988: if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsfunctionctxcontainer, ctx));
989: else *ctx = NULL;
990: }
991: PetscFunctionReturn(PETSC_SUCCESS);
992: }
994: /*@C
995: DMTSSetIJacobian - set `TS` Jacobian evaluation function into a `DMTS`
997: Not Collective
999: Input Parameters:
1000: + dm - `DM` to be used with `TS`
1001: . func - Jacobian evaluation routine, see `TSIJacobianFn` for the calling sequence
1002: - ctx - context for residual evaluation
1004: Level: developer
1006: Note:
1007: `TSSetIJacobian()` is normally used, but it calls this function internally because the user context is actually
1008: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
1009: not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1011: .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `TSIJacobianFn`, `DMTSGetIJacobian()`, `TSSetIJacobian()`
1012: @*/
1013: PetscErrorCode DMTSSetIJacobian(DM dm, TSIJacobianFn *func, void *ctx)
1014: {
1015: DMTS tsdm;
1017: PetscFunctionBegin;
1019: PetscCall(DMGetDMTSWrite(dm, &tsdm));
1020: if (func) tsdm->ops->ijacobian = func;
1021: if (ctx) {
1022: PetscContainer ctxcontainer;
1023: PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
1024: PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
1025: PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", (PetscObject)ctxcontainer));
1026: tsdm->ijacobianctxcontainer = ctxcontainer;
1027: PetscCall(PetscContainerDestroy(&ctxcontainer));
1028: }
1029: PetscFunctionReturn(PETSC_SUCCESS);
1030: }
1032: /*@C
1033: DMTSSetIJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function into a `DMTS`
1035: Not Collective
1037: Input Parameters:
1038: + dm - `DM` to be used with `TS`
1039: - f - Jacobian evaluation context destroy function
1041: Level: developer
1043: Note:
1044: `TSSetIJacobianContextDestroy()` is normally used, but it calls this function internally because the user context is actually
1045: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
1046: not.
1048: Developer Notes:
1049: If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1051: .seealso: [](ch_ts), `DMTS`, `TSSetIJacobianContextDestroy()`, `TSSetI2JacobianContextDestroy()`, `DMTSSetIJacobian()`, `TSSetIJacobian()`
1052: @*/
1053: PetscErrorCode DMTSSetIJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
1054: {
1055: DMTS tsdm;
1057: PetscFunctionBegin;
1059: PetscCall(DMGetDMTSWrite(dm, &tsdm));
1060: if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ijacobianctxcontainer, f));
1061: PetscFunctionReturn(PETSC_SUCCESS);
1062: }
1064: PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM dm)
1065: {
1066: DMTS tsdm;
1068: PetscFunctionBegin;
1070: PetscCall(DMGetDMTSWrite(dm, &tsdm));
1071: PetscCall(DMTSUnsetIJacobianContext_DMTS(tsdm));
1072: PetscFunctionReturn(PETSC_SUCCESS);
1073: }
1075: /*@C
1076: DMTSGetIJacobian - get `TS` Jacobian evaluation function from a `DMTS`
1078: Not Collective
1080: Input Parameter:
1081: . dm - `DM` to be used with `TS`
1083: Output Parameters:
1084: + func - Jacobian evaluation function, for calling sequence see `TSIJacobianFn`
1085: - ctx - context for residual evaluation
1087: Level: developer
1089: Note:
1090: `TSGetIJacobian()` is normally used, but it calls this function internally because the user context is actually
1091: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
1092: not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1094: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetIJacobian()`, `TSIJacobianFn`
1095: @*/
1096: PetscErrorCode DMTSGetIJacobian(DM dm, TSIJacobianFn **func, void **ctx)
1097: {
1098: DMTS tsdm;
1100: PetscFunctionBegin;
1102: PetscCall(DMGetDMTS(dm, &tsdm));
1103: if (func) *func = tsdm->ops->ijacobian;
1104: if (ctx) {
1105: if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ijacobianctxcontainer, ctx));
1106: else *ctx = NULL;
1107: }
1108: PetscFunctionReturn(PETSC_SUCCESS);
1109: }
1111: /*@C
1112: DMTSSetRHSJacobian - set `TS` Jacobian evaluation function into a `DMTS`
1114: Not Collective
1116: Input Parameters:
1117: + dm - `DM` to be used with `TS`
1118: . func - Jacobian evaluation routine, for calling sequence see `TSIJacobianFn`
1119: - ctx - context for residual evaluation
1121: Level: developer
1123: Note:
1124: `TSSetRHSJacobian()` is normally used, but it calls this function internally because the user context is actually
1125: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
1126: not.
1128: Developer Notes:
1129: If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1131: .seealso: [](ch_ts), `DMTS`, `TSRHSJacobianFn`, `DMTSGetRHSJacobian()`, `TSSetRHSJacobian()`
1132: @*/
1133: PetscErrorCode DMTSSetRHSJacobian(DM dm, TSRHSJacobianFn *func, void *ctx)
1134: {
1135: DMTS tsdm;
1137: PetscFunctionBegin;
1139: PetscCall(DMGetDMTSWrite(dm, &tsdm));
1140: if (func) tsdm->ops->rhsjacobian = func;
1141: if (ctx) {
1142: PetscContainer ctxcontainer;
1143: PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer));
1144: PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
1145: PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", (PetscObject)ctxcontainer));
1146: tsdm->rhsjacobianctxcontainer = ctxcontainer;
1147: PetscCall(PetscContainerDestroy(&ctxcontainer));
1148: }
1149: PetscFunctionReturn(PETSC_SUCCESS);
1150: }
1152: /*@C
1153: DMTSSetRHSJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function from a `DMTS`
1155: Not Collective
1157: Input Parameters:
1158: + dm - `DM` to be used with `TS`
1159: - f - Jacobian evaluation context destroy function
1161: Level: developer
1163: Note:
1164: The user usually calls `TSSetRHSJacobianContextDestroy()` which calls this routine
1166: .seealso: [](ch_ts), `DMTS`, `TS`, `TSSetRHSJacobianContextDestroy()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
1167: @*/
1168: PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
1169: {
1170: DMTS tsdm;
1172: PetscFunctionBegin;
1174: PetscCall(DMGetDMTSWrite(dm, &tsdm));
1175: if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsjacobianctxcontainer, f));
1176: PetscFunctionReturn(PETSC_SUCCESS);
1177: }
1179: PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM dm)
1180: {
1181: DMTS tsdm;
1183: PetscFunctionBegin;
1185: PetscCall(DMGetDMTSWrite(dm, &tsdm));
1186: PetscCall(DMTSUnsetRHSJacobianContext_DMTS(tsdm));
1187: PetscFunctionReturn(PETSC_SUCCESS);
1188: }
1190: /*@C
1191: DMTSGetRHSJacobian - get `TS` Jacobian evaluation function from a `DMTS`
1193: Not Collective
1195: Input Parameter:
1196: . dm - `DM` to be used with `TS`
1198: Output Parameters:
1199: + func - Jacobian evaluation function, for calling sequence see `TSRHSJacobianFn`
1200: - ctx - context for residual evaluation
1202: Level: developer
1204: Note:
1205: `TSGetRHSJacobian()` is normally used, but it calls this function internally because the user context is actually
1206: associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or
1207: not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.
1209: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetRHSJacobian()`, `TSRHSJacobianFn`
1210: @*/
1211: PetscErrorCode DMTSGetRHSJacobian(DM dm, TSRHSJacobianFn **func, void **ctx)
1212: {
1213: DMTS tsdm;
1215: PetscFunctionBegin;
1217: PetscCall(DMGetDMTS(dm, &tsdm));
1218: if (func) *func = tsdm->ops->rhsjacobian;
1219: if (ctx) {
1220: if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsjacobianctxcontainer, ctx));
1221: else *ctx = NULL;
1222: }
1223: PetscFunctionReturn(PETSC_SUCCESS);
1224: }
1226: /*@C
1227: DMTSSetIFunctionSerialize - sets functions used to view and load a `TSIFunctionFn` context
1229: Not Collective
1231: Input Parameters:
1232: + dm - `DM` to be used with `TS`
1233: . view - viewer function
1234: - load - loading function
1236: Level: developer
1238: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`
1239: @*/
1240: PetscErrorCode DMTSSetIFunctionSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1241: {
1242: DMTS tsdm;
1244: PetscFunctionBegin;
1246: PetscCall(DMGetDMTSWrite(dm, &tsdm));
1247: tsdm->ops->ifunctionview = view;
1248: tsdm->ops->ifunctionload = load;
1249: PetscFunctionReturn(PETSC_SUCCESS);
1250: }
1252: /*@C
1253: DMTSSetIJacobianSerialize - sets functions used to view and load a `TSIJacobianFn` context
1255: Not Collective
1257: Input Parameters:
1258: + dm - `DM` to be used with `TS`
1259: . view - viewer function
1260: - load - loading function
1262: Level: developer
1264: .seealso: [](ch_ts), `DMTS`, `DM`, `TS`
1265: @*/
1266: PetscErrorCode DMTSSetIJacobianSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1267: {
1268: DMTS tsdm;
1270: PetscFunctionBegin;
1272: PetscCall(DMGetDMTSWrite(dm, &tsdm));
1273: tsdm->ops->ijacobianview = view;
1274: tsdm->ops->ijacobianload = load;
1275: PetscFunctionReturn(PETSC_SUCCESS);
1276: }