1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
2 #include <petsc/private/dmimpl.h>
3
DMTSUnsetRHSFunctionContext_DMTS(DMTS tsdm)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 }
11
DMTSUnsetRHSJacobianContext_DMTS(DMTS tsdm)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 }
19
DMTSUnsetIFunctionContext_DMTS(DMTS tsdm)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 }
27
DMTSUnsetIJacobianContext_DMTS(DMTS tsdm)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 }
35
DMTSUnsetI2FunctionContext_DMTS(DMTS tsdm)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 }
43
DMTSUnsetI2JacobianContext_DMTS(DMTS tsdm)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 }
51
DMTSDestroy(DMTS * kdm)52 static PetscErrorCode DMTSDestroy(DMTS *kdm)
53 {
54 PetscFunctionBegin;
55 if (!*kdm) PetscFunctionReturn(PETSC_SUCCESS);
56 PetscValidHeaderSpecific(*kdm, DMTS_CLASSID, 1);
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 }
71
DMTSLoad(DMTS kdm,PetscViewer viewer)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 PetscCtx ctx;
80
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 PetscCtx ctx;
89
90 PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer, &ctx));
91 PetscCall((*kdm->ops->ijacobianload)(&ctx, viewer));
92 }
93 PetscFunctionReturn(PETSC_SUCCESS);
94 }
95
DMTSView(DMTS kdm,PetscViewer viewer)96 PetscErrorCode DMTSView(DMTS kdm, PetscViewer viewer)
97 {
98 PetscBool isascii, isbinary;
99
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;
106
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;
131
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 PetscCtx ctx;
140
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 PetscCtx ctx;
152
153 PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer, &ctx));
154 PetscCall((*kdm->ops->ijacobianview)(ctx, viewer));
155 }
156 }
157 PetscFunctionReturn(PETSC_SUCCESS);
158 }
159
DMTSCreate(MPI_Comm comm,DMTS * kdm)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 }
167
168 /* Attaches the DMTS to the coarse level.
169 * Under what conditions should we copy versus duplicate?
170 */
DMCoarsenHook_DMTS(DM dm,DM dmc,PetscCtx ctx)171 static PetscErrorCode DMCoarsenHook_DMTS(DM dm, DM dmc, PetscCtx ctx)
172 {
173 PetscFunctionBegin;
174 PetscCall(DMCopyDMTS(dm, dmc));
175 PetscFunctionReturn(PETSC_SUCCESS);
176 }
177
178 /* This could restrict auxiliary information to the coarse level.
179 */
DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,PetscCtx ctx)180 static PetscErrorCode DMRestrictHook_DMTS(DM dm, Mat Restrict, Vec rscale, Mat Inject, DM dmc, PetscCtx ctx)
181 {
182 PetscFunctionBegin;
183 PetscFunctionReturn(PETSC_SUCCESS);
184 }
185
DMSubDomainHook_DMTS(DM dm,DM subdm,PetscCtx ctx)186 static PetscErrorCode DMSubDomainHook_DMTS(DM dm, DM subdm, PetscCtx ctx)
187 {
188 PetscFunctionBegin;
189 PetscCall(DMCopyDMTS(dm, subdm));
190 PetscFunctionReturn(PETSC_SUCCESS);
191 }
192
193 /* This could restrict auxiliary information to the coarse level.
194 */
DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,PetscCtx ctx)195 static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, PetscCtx ctx)
196 {
197 PetscFunctionBegin;
198 PetscFunctionReturn(PETSC_SUCCESS);
199 }
200
201 /*@C
202 DMTSCopy - copies the information in a `DMTS` to another `DMTS`
203
204 Not Collective
205
206 Input Parameters:
207 + kdm - Original `DMTS`
208 - nkdm - `DMTS` to receive the data, should have been created with `DMTSCreate()`
209
210 Level: developer
211
212 .seealso: [](ch_ts), `DMTSCreate()`, `DMTSDestroy()`
213 @*/
DMTSCopy(DMTS kdm,DMTS nkdm)214 PetscErrorCode DMTSCopy(DMTS kdm, DMTS nkdm)
215 {
216 PetscFunctionBegin;
217 PetscValidHeaderSpecific(kdm, DMTS_CLASSID, 1);
218 PetscValidHeaderSpecific(nkdm, DMTS_CLASSID, 2);
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;
228
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));
242
243 nkdm->data = kdm->data;
244
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 */
250
251 /* implementation specific copy hooks */
252 PetscTryTypeMethod(kdm, duplicate, nkdm);
253 PetscFunctionReturn(PETSC_SUCCESS);
254 }
255
256 /*@C
257 DMGetDMTS - get read-only private `DMTS` context from a `DM`
258
259 Not Collective
260
261 Input Parameter:
262 . dm - `DM` to be used with `TS`
263
264 Output Parameter:
265 . tsdm - private `DMTS` context
266
267 Level: developer
268
269 Notes:
270 Use `DMGetDMTSWrite()` if write access is needed. The `DMTSSetXXX()` API should be used wherever possible.
271
272 .seealso: [](ch_ts), `DMTS`, `DMGetDMTSWrite()`
273 @*/
DMGetDMTS(DM dm,DMTS * tsdm)274 PetscErrorCode DMGetDMTS(DM dm, DMTS *tsdm)
275 {
276 PetscFunctionBegin;
277 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 }
289
290 /*@C
291 DMGetDMTSWrite - get write access to private `DMTS` context from a `DM`
292
293 Not Collective
294
295 Input Parameter:
296 . dm - `DM` to be used with `TS`
297
298 Output Parameter:
299 . tsdm - private `DMTS` context
300
301 Level: developer
302
303 .seealso: [](ch_ts), `DMTS`, `DMGetDMTS()`
304 @*/
DMGetDMTSWrite(DM dm,DMTS * tsdm)305 PetscErrorCode DMGetDMTSWrite(DM dm, DMTS *tsdm)
306 {
307 DMTS sdm;
308
309 PetscFunctionBegin;
310 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 }
325
326 /*@C
327 DMCopyDMTS - copies a `DMTS` context to a new `DM`
328
329 Logically Collective
330
331 Input Parameters:
332 + dmsrc - `DM` to obtain context from
333 - dmdest - `DM` to add context to
334
335 Level: developer
336
337 Note:
338 The context is copied by reference. This function does not ensure that a context exists.
339
340 .seealso: [](ch_ts), `DMTS`, `DMGetDMTS()`, `TSSetDM()`
341 @*/
DMCopyDMTS(DM dmsrc,DM dmdest)342 PetscErrorCode DMCopyDMTS(DM dmsrc, DM dmdest)
343 {
344 PetscFunctionBegin;
345 PetscValidHeaderSpecific(dmsrc, DM_CLASSID, 1);
346 PetscValidHeaderSpecific(dmdest, DM_CLASSID, 2);
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 }
354
355 /*@C
356 DMTSSetIFunction - set `TS` implicit function evaluation function into a `DMTS`
357
358 Not Collective
359
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
364
365 Level: developer
366
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.
371
372 .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `TSIFunctionFn`
373 @*/
DMTSSetIFunction(DM dm,TSIFunctionFn * func,PetscCtx ctx)374 PetscErrorCode DMTSSetIFunction(DM dm, TSIFunctionFn *func, PetscCtx ctx)
375 {
376 DMTS tsdm;
377
378 PetscFunctionBegin;
379 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 }
392
393 /*@C
394 DMTSSetIFunctionContextDestroy - set `TS` implicit evaluation context destroy function into a `DMTS`
395
396 Not Collective
397
398 Input Parameters:
399 + dm - `DM` to be used with `TS`
400 - f - implicit evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
401
402 Level: developer
403
404 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetIFunction()`, `TSSetIFunction()`, `PetscCtxDestroyFn`
405 @*/
DMTSSetIFunctionContextDestroy(DM dm,PetscCtxDestroyFn * f)406 PetscErrorCode DMTSSetIFunctionContextDestroy(DM dm, PetscCtxDestroyFn *f)
407 {
408 DMTS tsdm;
409
410 PetscFunctionBegin;
411 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
412 PetscCall(DMGetDMTSWrite(dm, &tsdm));
413 if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->ifunctionctxcontainer, f));
414 PetscFunctionReturn(PETSC_SUCCESS);
415 }
416
DMTSUnsetIFunctionContext_Internal(DM dm)417 PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM dm)
418 {
419 DMTS tsdm;
420
421 PetscFunctionBegin;
422 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
423 PetscCall(DMGetDMTSWrite(dm, &tsdm));
424 PetscCall(DMTSUnsetIFunctionContext_DMTS(tsdm));
425 PetscFunctionReturn(PETSC_SUCCESS);
426 }
427
428 /*@C
429 DMTSGetIFunction - get `TS` implicit residual evaluation function from a `DMTS`
430
431 Not Collective
432
433 Input Parameter:
434 . dm - `DM` to be used with `TS`
435
436 Output Parameters:
437 + func - function evaluation function, for calling sequence see `TSIFunctionFn`
438 - ctx - context for residual evaluation
439
440 Level: developer
441
442 Note:
443 `TSGetIFunction()` is normally used, but it calls this function internally because the user context is actually
444 associated with the `DM`.
445
446 .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `DMTSSetIFunction()`, `TSIFunctionFn`
447 @*/
DMTSGetIFunction(DM dm,TSIFunctionFn ** func,PetscCtxRt ctx)448 PetscErrorCode DMTSGetIFunction(DM dm, TSIFunctionFn **func, PetscCtxRt ctx)
449 {
450 DMTS tsdm;
451
452 PetscFunctionBegin;
453 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 *(void **)ctx = NULL;
459 }
460 PetscFunctionReturn(PETSC_SUCCESS);
461 }
462
463 /*@C
464 DMTSSetI2Function - set `TS` implicit function evaluation function for 2nd order systems into a `TSDM`
465
466 Not Collective
467
468 Input Parameters:
469 + dm - `DM` to be used with `TS`
470 . fun - function evaluation routine
471 - ctx - context for residual evaluation
472
473 Level: developer
474
475 Note:
476 `TSSetI2Function()` is normally used, but it calls this function internally because the user context is actually
477 associated with the `DM`.
478
479 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSSetI2Function()`
480 @*/
DMTSSetI2Function(DM dm,TSI2FunctionFn * fun,PetscCtx ctx)481 PetscErrorCode DMTSSetI2Function(DM dm, TSI2FunctionFn *fun, PetscCtx ctx)
482 {
483 DMTS tsdm;
484
485 PetscFunctionBegin;
486 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 }
499
500 /*@C
501 DMTSSetI2FunctionContextDestroy - set `TS` implicit evaluation for 2nd order systems context destroy into a `DMTS`
502
503 Not Collective
504
505 Input Parameters:
506 + dm - `DM` to be used with `TS`
507 - f - implicit evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
508
509 Level: developer
510
511 Note:
512 `TSSetI2FunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually
513 associated with the `DM`.
514
515 .seealso: [](ch_ts), `DMTS`, `TSSetI2FunctionContextDestroy()`, `DMTSSetI2Function()`, `TSSetI2Function()`
516 @*/
DMTSSetI2FunctionContextDestroy(DM dm,PetscCtxDestroyFn * f)517 PetscErrorCode DMTSSetI2FunctionContextDestroy(DM dm, PetscCtxDestroyFn *f)
518 {
519 DMTS tsdm;
520
521 PetscFunctionBegin;
522 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
523 PetscCall(DMGetDMTSWrite(dm, &tsdm));
524 if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->i2functionctxcontainer, f));
525 PetscFunctionReturn(PETSC_SUCCESS);
526 }
527
DMTSUnsetI2FunctionContext_Internal(DM dm)528 PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM dm)
529 {
530 DMTS tsdm;
531
532 PetscFunctionBegin;
533 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
534 PetscCall(DMGetDMTSWrite(dm, &tsdm));
535 PetscCall(DMTSUnsetI2FunctionContext_DMTS(tsdm));
536 PetscFunctionReturn(PETSC_SUCCESS);
537 }
538
539 /*@C
540 DMTSGetI2Function - get `TS` implicit residual evaluation function for 2nd order systems from a `DMTS`
541
542 Not Collective
543
544 Input Parameter:
545 . dm - `DM` to be used with `TS`
546
547 Output Parameters:
548 + fun - function evaluation function, for calling sequence see `TSSetI2Function()`
549 - ctx - context for residual evaluation
550
551 Level: developer
552
553 Note:
554 `TSGetI2Function()` is normally used, but it calls this function internally because the user context is actually
555 associated with the `DM`.
556
557 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetI2Function()`, `TSGetI2Function()`
558 @*/
DMTSGetI2Function(DM dm,TSI2FunctionFn ** fun,PetscCtxRt ctx)559 PetscErrorCode DMTSGetI2Function(DM dm, TSI2FunctionFn **fun, PetscCtxRt ctx)
560 {
561 DMTS tsdm;
562
563 PetscFunctionBegin;
564 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 *(void **)ctx = NULL;
570 }
571 PetscFunctionReturn(PETSC_SUCCESS);
572 }
573
574 /*@C
575 DMTSSetI2Jacobian - set `TS` implicit Jacobian evaluation function for 2nd order systems from a `DMTS`
576
577 Not Collective
578
579 Input Parameters:
580 + dm - `DM` to be used with `TS`
581 . jac - Jacobian evaluation routine
582 - ctx - context for Jacobian evaluation
583
584 Level: developer
585
586 Note:
587 `TSSetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
588 associated with the `DM`.
589
590 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSI2JacobianFn`, `TSSetI2Jacobian()`
591 @*/
DMTSSetI2Jacobian(DM dm,TSI2JacobianFn * jac,PetscCtx ctx)592 PetscErrorCode DMTSSetI2Jacobian(DM dm, TSI2JacobianFn *jac, PetscCtx ctx)
593 {
594 DMTS tsdm;
595
596 PetscFunctionBegin;
597 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 }
610
611 /*@C
612 DMTSSetI2JacobianContextDestroy - set `TS` implicit Jacobian evaluation for 2nd order systems context destroy function into a `DMTS`
613
614 Not Collective
615
616 Input Parameters:
617 + dm - `DM` to be used with `TS`
618 - f - implicit Jacobian evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
619
620 Level: developer
621
622 Note:
623 Normally `TSSetI2JacobianContextDestroy()` is used
624
625 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSSetI2JacobianContextDestroy()`, `DMTSSetI2Jacobian()`, `TSSetI2Jacobian()`
626 @*/
DMTSSetI2JacobianContextDestroy(DM dm,PetscCtxDestroyFn * f)627 PetscErrorCode DMTSSetI2JacobianContextDestroy(DM dm, PetscCtxDestroyFn *f)
628 {
629 DMTS tsdm;
630
631 PetscFunctionBegin;
632 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
633 PetscCall(DMGetDMTSWrite(dm, &tsdm));
634 if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->i2jacobianctxcontainer, f));
635 PetscFunctionReturn(PETSC_SUCCESS);
636 }
637
DMTSUnsetI2JacobianContext_Internal(DM dm)638 PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM dm)
639 {
640 DMTS tsdm;
641
642 PetscFunctionBegin;
643 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
644 PetscCall(DMGetDMTSWrite(dm, &tsdm));
645 PetscCall(DMTSUnsetI2JacobianContext_DMTS(tsdm));
646 PetscFunctionReturn(PETSC_SUCCESS);
647 }
648
649 /*@C
650 DMTSGetI2Jacobian - get `TS` implicit Jacobian evaluation function for 2nd order systems from a `DMTS`
651
652 Not Collective
653
654 Input Parameter:
655 . dm - `DM` to be used with `TS`
656
657 Output Parameters:
658 + jac - Jacobian evaluation function, for calling sequence see `TSI2JacobianFn`
659 - ctx - context for Jacobian evaluation
660
661 Level: developer
662
663 Note:
664 `TSGetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually
665 associated with the `DM`.
666
667 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()`, `TSI2JacobianFn`
668 @*/
DMTSGetI2Jacobian(DM dm,TSI2JacobianFn ** jac,PetscCtxRt ctx)669 PetscErrorCode DMTSGetI2Jacobian(DM dm, TSI2JacobianFn **jac, PetscCtxRt ctx)
670 {
671 DMTS tsdm;
672
673 PetscFunctionBegin;
674 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 *(void **)ctx = NULL;
680 }
681 PetscFunctionReturn(PETSC_SUCCESS);
682 }
683
684 /*@C
685 DMTSSetRHSFunction - set `TS` explicit residual evaluation function into a `DMTS`
686
687 Not Collective
688
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
693
694 Level: developer
695
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.
700
701 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSRHSFunctionFn`
702 @*/
DMTSSetRHSFunction(DM dm,TSRHSFunctionFn * func,PetscCtx ctx)703 PetscErrorCode DMTSSetRHSFunction(DM dm, TSRHSFunctionFn *func, PetscCtx ctx)
704 {
705 DMTS tsdm;
706
707 PetscFunctionBegin;
708 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 }
721
722 /*@C
723 DMTSSetRHSFunctionContextDestroy - set `TS` explicit residual evaluation context destroy function into a `DMTS`
724
725 Not Collective
726
727 Input Parameters:
728 + dm - `DM` to be used with `TS`
729 - f - explicit evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
730
731 Level: developer
732
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.
737
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.
740
741 .seealso: [](ch_ts), `DMTS`, `TSSetRHSFunctionContextDestroy()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()`
742 @*/
DMTSSetRHSFunctionContextDestroy(DM dm,PetscCtxDestroyFn * f)743 PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM dm, PetscCtxDestroyFn *f)
744 {
745 DMTS tsdm;
746
747 PetscFunctionBegin;
748 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
749 PetscCall(DMGetDMTSWrite(dm, &tsdm));
750 if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->rhsfunctionctxcontainer, f));
751 PetscFunctionReturn(PETSC_SUCCESS);
752 }
753
DMTSUnsetRHSFunctionContext_Internal(DM dm)754 PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM dm)
755 {
756 DMTS tsdm;
757
758 PetscFunctionBegin;
759 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
760 PetscCall(DMGetDMTSWrite(dm, &tsdm));
761 PetscCall(DMTSUnsetRHSFunctionContext_DMTS(tsdm));
762 tsdm->rhsfunctionctxcontainer = NULL;
763 PetscFunctionReturn(PETSC_SUCCESS);
764 }
765
766 /*@C
767 DMTSSetTransientVariable - sets function to transform from state to transient variables into a `DMTS`
768
769 Logically Collective
770
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
775
776 Level: developer
777
778 Notes:
779 Normally `TSSetTransientVariable()` is used
780
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
786
787 $$
788 dF/dP + shift * dF/dCdot dC/dP.
789 $$
790
791 .seealso: [](ch_ts), `DMTS`, `TS`, `TSBDF`, `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()`, `TSTransientVariableFn`
792 @*/
DMTSSetTransientVariable(DM dm,TSTransientVariableFn * tvar,PetscCtx ctx)793 PetscErrorCode DMTSSetTransientVariable(DM dm, TSTransientVariableFn *tvar, PetscCtx ctx)
794 {
795 DMTS dmts;
796
797 PetscFunctionBegin;
798 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
799 PetscCall(DMGetDMTSWrite(dm, &dmts));
800 dmts->ops->transientvar = tvar;
801 dmts->transientvarctx = ctx;
802 PetscFunctionReturn(PETSC_SUCCESS);
803 }
804
805 /*@C
806 DMTSGetTransientVariable - gets function to transform from state to transient variables set with `DMTSSetTransientVariable()` from a `TSDM`
807
808 Logically Collective
809
810 Input Parameter:
811 . dm - `DM` to be used with `TS`
812
813 Output Parameters:
814 + tvar - a function that transforms to transient variables, see `TSTransientVariableFn` for the calling sequence
815 - ctx - a context for tvar
816
817 Level: developer
818
819 Note:
820 Normally `TSSetTransientVariable()` is used
821
822 .seealso: [](ch_ts), `DMTS`, `DM`, `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()`, `TSTransientVariableFn`
823 @*/
DMTSGetTransientVariable(DM dm,TSTransientVariableFn ** tvar,PetscCtx ctx)824 PetscErrorCode DMTSGetTransientVariable(DM dm, TSTransientVariableFn **tvar, PetscCtx ctx)
825 {
826 DMTS dmts;
827
828 PetscFunctionBegin;
829 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
830 PetscCall(DMGetDMTS(dm, &dmts));
831 if (tvar) *tvar = dmts->ops->transientvar;
832 if (ctx) *(void **)ctx = dmts->transientvarctx;
833 PetscFunctionReturn(PETSC_SUCCESS);
834 }
835
836 /*@C
837 DMTSGetSolutionFunction - gets the `TS` solution evaluation function from a `DMTS`
838
839 Not Collective
840
841 Input Parameter:
842 . dm - `DM` to be used with `TS`
843
844 Output Parameters:
845 + func - solution function evaluation function, for calling sequence see `TSSolutionFn`
846 - ctx - context for solution evaluation
847
848 Level: developer
849
850 .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `DMTSSetSolutionFunction()`, `TSSolutionFn`
851 @*/
DMTSGetSolutionFunction(DM dm,TSSolutionFn ** func,PetscCtxRt ctx)852 PetscErrorCode DMTSGetSolutionFunction(DM dm, TSSolutionFn **func, PetscCtxRt ctx)
853 {
854 DMTS tsdm;
855
856 PetscFunctionBegin;
857 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
858 PetscCall(DMGetDMTS(dm, &tsdm));
859 if (func) *func = tsdm->ops->solution;
860 if (ctx) *(void **)ctx = tsdm->solutionctx;
861 PetscFunctionReturn(PETSC_SUCCESS);
862 }
863
864 /*@C
865 DMTSSetSolutionFunction - set `TS` solution evaluation function into a `DMTS`
866
867 Not Collective
868
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
873
874 Level: developer
875
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.
880
881 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSGetSolutionFunction()`, `TSSolutionFn`
882 @*/
DMTSSetSolutionFunction(DM dm,TSSolutionFn * func,PetscCtx ctx)883 PetscErrorCode DMTSSetSolutionFunction(DM dm, TSSolutionFn *func, PetscCtx ctx)
884 {
885 DMTS tsdm;
886
887 PetscFunctionBegin;
888 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
889 PetscCall(DMGetDMTSWrite(dm, &tsdm));
890 if (func) tsdm->ops->solution = func;
891 if (ctx) tsdm->solutionctx = ctx;
892 PetscFunctionReturn(PETSC_SUCCESS);
893 }
894
895 /*@C
896 DMTSSetForcingFunction - set `TS` forcing function evaluation function into a `DMTS`
897
898 Not Collective
899
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
904
905 Level: developer
906
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.
911
912 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSForcingFn`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()`
913 @*/
DMTSSetForcingFunction(DM dm,TSForcingFn * func,PetscCtx ctx)914 PetscErrorCode DMTSSetForcingFunction(DM dm, TSForcingFn *func, PetscCtx ctx)
915 {
916 DMTS tsdm;
917
918 PetscFunctionBegin;
919 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
920 PetscCall(DMGetDMTSWrite(dm, &tsdm));
921 if (func) tsdm->ops->forcing = func;
922 if (ctx) tsdm->forcingctx = ctx;
923 PetscFunctionReturn(PETSC_SUCCESS);
924 }
925
926 /*@C
927 DMTSGetForcingFunction - get `TS` forcing function evaluation function from a `DMTS`
928
929 Not Collective
930
931 Input Parameter:
932 . dm - `DM` to be used with `TS`
933
934 Output Parameters:
935 + f - forcing function evaluation function; see `TSForcingFn` for the calling sequence
936 - ctx - context for solution evaluation
937
938 Level: developer
939
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.
944
945 .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `TSSetForcingFunction()`, `TSForcingFn`
946 @*/
DMTSGetForcingFunction(DM dm,TSForcingFn ** f,PetscCtxRt ctx)947 PetscErrorCode DMTSGetForcingFunction(DM dm, TSForcingFn **f, PetscCtxRt ctx)
948 {
949 DMTS tsdm;
950
951 PetscFunctionBegin;
952 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
953 PetscCall(DMGetDMTSWrite(dm, &tsdm));
954 if (f) *f = tsdm->ops->forcing;
955 if (ctx) *(void **)ctx = tsdm->forcingctx;
956 PetscFunctionReturn(PETSC_SUCCESS);
957 }
958
959 /*@C
960 DMTSGetRHSFunction - get `TS` explicit residual evaluation function from a `DMTS`
961
962 Not Collective
963
964 Input Parameter:
965 . dm - `DM` to be used with `TS`
966
967 Output Parameters:
968 + func - residual evaluation function, for calling sequence see `TSRHSFunctionFn`
969 - ctx - context for residual evaluation
970
971 Level: developer
972
973 Note:
974 `TSGetRHSFunction()` is normally used, but it calls this function internally because the user context is actually
975 associated with the DM.
976
977 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `TSRHSFunctionFn`, `TSGetRHSFunction()`
978 @*/
DMTSGetRHSFunction(DM dm,TSRHSFunctionFn ** func,PetscCtxRt ctx)979 PetscErrorCode DMTSGetRHSFunction(DM dm, TSRHSFunctionFn **func, PetscCtxRt ctx)
980 {
981 DMTS tsdm;
982
983 PetscFunctionBegin;
984 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 *(void **)ctx = NULL;
990 }
991 PetscFunctionReturn(PETSC_SUCCESS);
992 }
993
994 /*@C
995 DMTSSetIJacobian - set `TS` Jacobian evaluation function into a `DMTS`
996
997 Not Collective
998
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
1003
1004 Level: developer
1005
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.
1010
1011 .seealso: [](ch_ts), `DMTS`, `TS`, `DM`, `TSIJacobianFn`, `DMTSGetIJacobian()`, `TSSetIJacobian()`
1012 @*/
DMTSSetIJacobian(DM dm,TSIJacobianFn * func,PetscCtx ctx)1013 PetscErrorCode DMTSSetIJacobian(DM dm, TSIJacobianFn *func, PetscCtx ctx)
1014 {
1015 DMTS tsdm;
1016
1017 PetscFunctionBegin;
1018 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 }
1031
1032 /*@C
1033 DMTSSetIJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function into a `DMTS`
1034
1035 Not Collective
1036
1037 Input Parameters:
1038 + dm - `DM` to be used with `TS`
1039 - f - Jacobian evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
1040
1041 Level: developer
1042
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.
1047
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.
1050
1051 .seealso: [](ch_ts), `DMTS`, `TSSetIJacobianContextDestroy()`, `TSSetI2JacobianContextDestroy()`, `DMTSSetIJacobian()`, `TSSetIJacobian()`
1052 @*/
DMTSSetIJacobianContextDestroy(DM dm,PetscCtxDestroyFn * f)1053 PetscErrorCode DMTSSetIJacobianContextDestroy(DM dm, PetscCtxDestroyFn *f)
1054 {
1055 DMTS tsdm;
1056
1057 PetscFunctionBegin;
1058 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1059 PetscCall(DMGetDMTSWrite(dm, &tsdm));
1060 if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->ijacobianctxcontainer, f));
1061 PetscFunctionReturn(PETSC_SUCCESS);
1062 }
1063
DMTSUnsetIJacobianContext_Internal(DM dm)1064 PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM dm)
1065 {
1066 DMTS tsdm;
1067
1068 PetscFunctionBegin;
1069 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1070 PetscCall(DMGetDMTSWrite(dm, &tsdm));
1071 PetscCall(DMTSUnsetIJacobianContext_DMTS(tsdm));
1072 PetscFunctionReturn(PETSC_SUCCESS);
1073 }
1074
1075 /*@C
1076 DMTSGetIJacobian - get `TS` Jacobian evaluation function from a `DMTS`
1077
1078 Not Collective
1079
1080 Input Parameter:
1081 . dm - `DM` to be used with `TS`
1082
1083 Output Parameters:
1084 + func - Jacobian evaluation function, for calling sequence see `TSIJacobianFn`
1085 - ctx - context for residual evaluation
1086
1087 Level: developer
1088
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.
1093
1094 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetIJacobian()`, `TSIJacobianFn`
1095 @*/
DMTSGetIJacobian(DM dm,TSIJacobianFn ** func,PetscCtxRt ctx)1096 PetscErrorCode DMTSGetIJacobian(DM dm, TSIJacobianFn **func, PetscCtxRt ctx)
1097 {
1098 DMTS tsdm;
1099
1100 PetscFunctionBegin;
1101 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 *(void **)ctx = NULL;
1107 }
1108 PetscFunctionReturn(PETSC_SUCCESS);
1109 }
1110
1111 /*@C
1112 DMTSSetRHSJacobian - set `TS` Jacobian evaluation function into a `DMTS`
1113
1114 Not Collective
1115
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
1120
1121 Level: developer
1122
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.
1127
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.
1130
1131 .seealso: [](ch_ts), `DMTS`, `TSRHSJacobianFn`, `DMTSGetRHSJacobian()`, `TSSetRHSJacobian()`
1132 @*/
DMTSSetRHSJacobian(DM dm,TSRHSJacobianFn * func,PetscCtx ctx)1133 PetscErrorCode DMTSSetRHSJacobian(DM dm, TSRHSJacobianFn *func, PetscCtx ctx)
1134 {
1135 DMTS tsdm;
1136
1137 PetscFunctionBegin;
1138 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 }
1151
1152 /*@C
1153 DMTSSetRHSJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function from a `DMTS`
1154
1155 Not Collective
1156
1157 Input Parameters:
1158 + dm - `DM` to be used with `TS`
1159 - f - Jacobian evaluation context destroy function, see `PetscCtxDestroyFn` for its calling sequence
1160
1161 Level: developer
1162
1163 Note:
1164 The user usually calls `TSSetRHSJacobianContextDestroy()` which calls this routine
1165
1166 .seealso: [](ch_ts), `DMTS`, `TS`, `TSSetRHSJacobianContextDestroy()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()`
1167 @*/
DMTSSetRHSJacobianContextDestroy(DM dm,PetscCtxDestroyFn * f)1168 PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM dm, PetscCtxDestroyFn *f)
1169 {
1170 DMTS tsdm;
1171
1172 PetscFunctionBegin;
1173 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1174 PetscCall(DMGetDMTSWrite(dm, &tsdm));
1175 if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerSetCtxDestroy(tsdm->rhsjacobianctxcontainer, f));
1176 PetscFunctionReturn(PETSC_SUCCESS);
1177 }
1178
DMTSUnsetRHSJacobianContext_Internal(DM dm)1179 PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM dm)
1180 {
1181 DMTS tsdm;
1182
1183 PetscFunctionBegin;
1184 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1185 PetscCall(DMGetDMTSWrite(dm, &tsdm));
1186 PetscCall(DMTSUnsetRHSJacobianContext_DMTS(tsdm));
1187 PetscFunctionReturn(PETSC_SUCCESS);
1188 }
1189
1190 /*@C
1191 DMTSGetRHSJacobian - get `TS` Jacobian evaluation function from a `DMTS`
1192
1193 Not Collective
1194
1195 Input Parameter:
1196 . dm - `DM` to be used with `TS`
1197
1198 Output Parameters:
1199 + func - Jacobian evaluation function, for calling sequence see `TSRHSJacobianFn`
1200 - ctx - context for residual evaluation
1201
1202 Level: developer
1203
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.
1208
1209 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`, `DMTSSetRHSJacobian()`, `TSRHSJacobianFn`
1210 @*/
DMTSGetRHSJacobian(DM dm,TSRHSJacobianFn ** func,PetscCtxRt ctx)1211 PetscErrorCode DMTSGetRHSJacobian(DM dm, TSRHSJacobianFn **func, PetscCtxRt ctx)
1212 {
1213 DMTS tsdm;
1214
1215 PetscFunctionBegin;
1216 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
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 *(void **)ctx = NULL;
1222 }
1223 PetscFunctionReturn(PETSC_SUCCESS);
1224 }
1225
1226 /*@C
1227 DMTSSetIFunctionSerialize - sets functions used to view and load a `TSIFunctionFn` context
1228
1229 Not Collective
1230
1231 Input Parameters:
1232 + dm - `DM` to be used with `TS`
1233 . view - viewer function
1234 - load - loading function
1235
1236 Level: developer
1237
1238 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`
1239 @*/
DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (* view)(void *,PetscViewer),PetscErrorCode (* load)(void **,PetscViewer))1240 PetscErrorCode DMTSSetIFunctionSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1241 {
1242 DMTS tsdm;
1243
1244 PetscFunctionBegin;
1245 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1246 PetscCall(DMGetDMTSWrite(dm, &tsdm));
1247 tsdm->ops->ifunctionview = view;
1248 tsdm->ops->ifunctionload = load;
1249 PetscFunctionReturn(PETSC_SUCCESS);
1250 }
1251
1252 /*@C
1253 DMTSSetIJacobianSerialize - sets functions used to view and load a `TSIJacobianFn` context
1254
1255 Not Collective
1256
1257 Input Parameters:
1258 + dm - `DM` to be used with `TS`
1259 . view - viewer function
1260 - load - loading function
1261
1262 Level: developer
1263
1264 .seealso: [](ch_ts), `DMTS`, `DM`, `TS`
1265 @*/
DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (* view)(void *,PetscViewer),PetscErrorCode (* load)(void **,PetscViewer))1266 PetscErrorCode DMTSSetIJacobianSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer))
1267 {
1268 DMTS tsdm;
1269
1270 PetscFunctionBegin;
1271 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1272 PetscCall(DMGetDMTSWrite(dm, &tsdm));
1273 tsdm->ops->ijacobianview = view;
1274 tsdm->ops->ijacobianload = load;
1275 PetscFunctionReturn(PETSC_SUCCESS);
1276 }
1277