xref: /petsc/src/ts/utils/dmts.c (revision 9f537349bd91d7333cf178eb5d8f8d6f36dba08d)
1 #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
2 #include <petsc/private/dmimpl.h>
3 
4 static PetscErrorCode DMTSDestroy(DMTS *kdm)
5 {
6   PetscErrorCode ierr;
7 
8   PetscFunctionBegin;
9   if (!*kdm) PetscFunctionReturn(0);
10   PetscValidHeaderSpecific((*kdm),DMTS_CLASSID,1);
11   if (--((PetscObject)(*kdm))->refct > 0) {*kdm = 0; PetscFunctionReturn(0);}
12   if ((*kdm)->ops->destroy) {ierr = ((*kdm)->ops->destroy)(*kdm);CHKERRQ(ierr);}
13   ierr = PetscHeaderDestroy(kdm);CHKERRQ(ierr);
14   PetscFunctionReturn(0);
15 }
16 
17 PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer)
18 {
19   PetscErrorCode ierr;
20 
21   PetscFunctionBegin;
22   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
23   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
24   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
25   if (kdm->ops->ifunctionload) {
26     ierr = (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);CHKERRQ(ierr);
27   }
28   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
29   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
30   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr);
31   if (kdm->ops->ijacobianload) {
32     ierr = (*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer);CHKERRQ(ierr);
33   }
34   PetscFunctionReturn(0);
35 }
36 
37 PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
38 {
39   PetscErrorCode ierr;
40   PetscBool      isascii,isbinary;
41 
42   PetscFunctionBegin;
43   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
44   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
45   if (isascii) {
46 #if defined(PETSC_SERIALIZE_FUNCTIONS)
47     const char *fname;
48 
49     ierr = PetscFPTFind(kdm->ops->ifunction,&fname);CHKERRQ(ierr);
50     if (fname) {
51       ierr = PetscViewerASCIIPrintf(viewer,"  IFunction used by TS: %s\n",fname);CHKERRQ(ierr);
52     }
53     ierr = PetscFPTFind(kdm->ops->ijacobian,&fname);CHKERRQ(ierr);
54     if (fname) {
55       ierr = PetscViewerASCIIPrintf(viewer,"  IJacobian function used by TS: %s\n",fname);CHKERRQ(ierr);
56     }
57 #endif
58   } else if (isbinary) {
59     struct {
60       TSIFunction ifunction;
61     } funcstruct;
62     struct {
63       PetscErrorCode (*ifunctionview)(void*,PetscViewer);
64     } funcviewstruct;
65     struct {
66       PetscErrorCode (*ifunctionload)(void**,PetscViewer);
67     } funcloadstruct;
68     struct {
69       TSIJacobian ijacobian;
70     } jacstruct;
71     struct {
72       PetscErrorCode (*ijacobianview)(void*,PetscViewer);
73     } jacviewstruct;
74     struct {
75       PetscErrorCode (*ijacobianload)(void**,PetscViewer);
76     } jacloadstruct;
77 
78     funcstruct.ifunction         = kdm->ops->ifunction;
79     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
80     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
81     ierr = PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
82     ierr = PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
83     ierr = PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
84     if (kdm->ops->ifunctionview) {
85       ierr = (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);CHKERRQ(ierr);
86     }
87     jacstruct.ijacobian = kdm->ops->ijacobian;
88     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
89     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
90     ierr = PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
91     ierr = PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
92     ierr = PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
93     if (kdm->ops->ijacobianview) {
94       ierr = (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);CHKERRQ(ierr);
95     }
96   }
97   PetscFunctionReturn(0);
98 }
99 
100 static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
101 {
102   PetscErrorCode ierr;
103 
104   PetscFunctionBegin;
105   ierr = TSInitializePackage();CHKERRQ(ierr);
106   ierr = PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);CHKERRQ(ierr);
107   PetscFunctionReturn(0);
108 }
109 
110 /* Attaches the DMTS to the coarse level.
111  * Under what conditions should we copy versus duplicate?
112  */
113 static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
114 {
115   PetscErrorCode ierr;
116 
117   PetscFunctionBegin;
118   ierr = DMCopyDMTS(dm,dmc);CHKERRQ(ierr);
119   PetscFunctionReturn(0);
120 }
121 
122 /* This could restrict auxiliary information to the coarse level.
123  */
124 static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
125 {
126 
127   PetscFunctionBegin;
128   PetscFunctionReturn(0);
129 }
130 
131 static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
132 {
133   PetscErrorCode ierr;
134 
135   PetscFunctionBegin;
136   ierr = DMCopyDMTS(dm,subdm);CHKERRQ(ierr);
137   PetscFunctionReturn(0);
138 }
139 
140 /* This could restrict auxiliary information to the coarse level.
141  */
142 static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
143 {
144   PetscFunctionBegin;
145   PetscFunctionReturn(0);
146 }
147 
148 /*@C
149    DMTSCopy - copies the information in a DMTS to another DMTS
150 
151    Not Collective
152 
153    Input Argument:
154 +  kdm - Original DMTS
155 -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()
156 
157    Level: developer
158 
159 .seealso: DMTSCreate(), DMTSDestroy()
160 @*/
161 PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
162 {
163   PetscErrorCode ierr;
164 
165   PetscFunctionBegin;
166   PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1);
167   PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2);
168   nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
169   nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
170   nkdm->ops->ifunction   = kdm->ops->ifunction;
171   nkdm->ops->ijacobian   = kdm->ops->ijacobian;
172   nkdm->ops->i2function  = kdm->ops->i2function;
173   nkdm->ops->i2jacobian  = kdm->ops->i2jacobian;
174   nkdm->ops->solution    = kdm->ops->solution;
175   nkdm->ops->destroy     = kdm->ops->destroy;
176   nkdm->ops->duplicate   = kdm->ops->duplicate;
177 
178   nkdm->rhsfunctionctx = kdm->rhsfunctionctx;
179   nkdm->rhsjacobianctx = kdm->rhsjacobianctx;
180   nkdm->ifunctionctx   = kdm->ifunctionctx;
181   nkdm->ijacobianctx   = kdm->ijacobianctx;
182   nkdm->i2functionctx  = kdm->i2functionctx;
183   nkdm->i2jacobianctx  = kdm->i2jacobianctx;
184   nkdm->solutionctx    = kdm->solutionctx;
185 
186   nkdm->data = kdm->data;
187 
188   /*
189   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
190   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
191   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
192   */
193 
194   /* implementation specific copy hooks */
195   if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);}
196   PetscFunctionReturn(0);
197 }
198 
199 /*@C
200    DMGetDMTS - get read-only private DMTS context from a DM
201 
202    Not Collective
203 
204    Input Argument:
205 .  dm - DM to be used with TS
206 
207    Output Argument:
208 .  tsdm - private DMTS context
209 
210    Level: developer
211 
212    Notes:
213    Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
214 
215 .seealso: DMGetDMTSWrite()
216 @*/
217 PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
218 {
219   PetscErrorCode ierr;
220 
221   PetscFunctionBegin;
222   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
223   *tsdm = (DMTS) dm->dmts;
224   if (!*tsdm) {
225     ierr = PetscInfo(dm,"Creating new DMTS\n");CHKERRQ(ierr);
226     ierr = DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm);CHKERRQ(ierr);
227     dm->dmts = (PetscObject) *tsdm;
228     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr);
229     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr);
230   }
231   PetscFunctionReturn(0);
232 }
233 
234 /*@C
235    DMGetDMTSWrite - get write access to private DMTS context from a DM
236 
237    Not Collective
238 
239    Input Argument:
240 .  dm - DM to be used with TS
241 
242    Output Argument:
243 .  tsdm - private DMTS context
244 
245    Level: developer
246 
247 .seealso: DMGetDMTS()
248 @*/
249 PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
250 {
251   PetscErrorCode ierr;
252   DMTS           sdm;
253 
254   PetscFunctionBegin;
255   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
256   ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr);
257   if (!sdm->originaldm) sdm->originaldm = dm;
258   if (sdm->originaldm != dm) {  /* Copy on write */
259     DMTS oldsdm = sdm;
260     ierr     = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr);
261     ierr     = DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);CHKERRQ(ierr);
262     ierr     = DMTSCopy(oldsdm,sdm);CHKERRQ(ierr);
263     ierr     = DMTSDestroy((DMTS*)&dm->dmts);CHKERRQ(ierr);
264     dm->dmts = (PetscObject) sdm;
265   }
266   *tsdm = sdm;
267   PetscFunctionReturn(0);
268 }
269 
270 /*@C
271    DMCopyDMTS - copies a DM context to a new DM
272 
273    Logically Collective
274 
275    Input Arguments:
276 +  dmsrc - DM to obtain context from
277 -  dmdest - DM to add context to
278 
279    Level: developer
280 
281    Note:
282    The context is copied by reference. This function does not ensure that a context exists.
283 
284 .seealso: DMGetDMTS(), TSSetDM()
285 @*/
286 PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
287 {
288   PetscErrorCode ierr;
289 
290   PetscFunctionBegin;
291   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
292   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
293   ierr         = DMTSDestroy((DMTS*)&dmdest->dmts);CHKERRQ(ierr);
294   dmdest->dmts = dmsrc->dmts;
295   ierr         = PetscObjectReference(dmdest->dmts);CHKERRQ(ierr);
296   ierr         = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr);
297   ierr         = DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr);
298   PetscFunctionReturn(0);
299 }
300 
301 /*@C
302    DMTSSetIFunction - set TS implicit function evaluation function
303 
304    Not Collective
305 
306    Input Arguments:
307 +  dm - DM to be used with TS
308 .  func - function evaluation function, see TSSetIFunction() for calling sequence
309 -  ctx - context for residual evaluation
310 
311    Level: advanced
312 
313    Note:
314    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
315    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
316    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
317 
318 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
319 @*/
320 PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
321 {
322   PetscErrorCode ierr;
323   DMTS           tsdm;
324 
325   PetscFunctionBegin;
326   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
327   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
328   if (func) tsdm->ops->ifunction = func;
329   if (ctx)  tsdm->ifunctionctx = ctx;
330   PetscFunctionReturn(0);
331 }
332 
333 /*@C
334    DMTSGetIFunction - get TS implicit residual evaluation function
335 
336    Not Collective
337 
338    Input Argument:
339 .  dm - DM to be used with TS
340 
341    Output Arguments:
342 +  func - function evaluation function, see TSSetIFunction() for calling sequence
343 -  ctx - context for residual evaluation
344 
345    Level: advanced
346 
347    Note:
348    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
349    associated with the DM.
350 
351 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
352 @*/
353 PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
354 {
355   PetscErrorCode ierr;
356   DMTS           tsdm;
357 
358   PetscFunctionBegin;
359   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
360   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
361   if (func) *func = tsdm->ops->ifunction;
362   if (ctx)  *ctx = tsdm->ifunctionctx;
363   PetscFunctionReturn(0);
364 }
365 
366 /*@C
367    DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems
368 
369    Not Collective
370 
371    Input Arguments:
372 +  dm - DM to be used with TS
373 .  fun - function evaluation function, see TSSetI2Function() for calling sequence
374 -  ctx - context for residual evaluation
375 
376    Level: advanced
377 
378    Note:
379    TSSetI2Function() is normally used, but it calls this function internally because the user context is actually
380    associated with the DM.
381 
382 .seealso: TSSetI2Function()
383 @*/
384 PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx)
385 {
386   DMTS           tsdm;
387   PetscErrorCode ierr;
388 
389   PetscFunctionBegin;
390   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
391   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
392   if (fun) tsdm->ops->i2function = fun;
393   if (ctx) tsdm->i2functionctx   = ctx;
394   PetscFunctionReturn(0);
395 }
396 
397 /*@C
398    DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems
399 
400    Not Collective
401 
402    Input Argument:
403 .  dm - DM to be used with TS
404 
405    Output Arguments:
406 +  fun - function evaluation function, see TSSetI2Function() for calling sequence
407 -  ctx - context for residual evaluation
408 
409    Level: advanced
410 
411    Note:
412    TSGetI2Function() is normally used, but it calls this function internally because the user context is actually
413    associated with the DM.
414 
415 .seealso: DMTSSetI2Function(),TSGetI2Function()
416 @*/
417 PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx)
418 {
419   DMTS           tsdm;
420   PetscErrorCode ierr;
421 
422   PetscFunctionBegin;
423   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
424   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
425   if (fun) *fun = tsdm->ops->i2function;
426   if (ctx) *ctx = tsdm->i2functionctx;
427   PetscFunctionReturn(0);
428 }
429 
430 /*@C
431    DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems
432 
433    Not Collective
434 
435    Input Arguments:
436 +  dm - DM to be used with TS
437 .  fun - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
438 -  ctx - context for Jacobian evaluation
439 
440    Level: advanced
441 
442    Note:
443    TSSetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
444    associated with the DM.
445 
446 .seealso: TSSetI2Jacobian()
447 @*/
448 PetscErrorCode DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void *ctx)
449 {
450   DMTS           tsdm;
451   PetscErrorCode ierr;
452 
453   PetscFunctionBegin;
454   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
455   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
456   if (jac) tsdm->ops->i2jacobian = jac;
457   if (ctx) tsdm->i2jacobianctx   = ctx;
458   PetscFunctionReturn(0);
459 }
460 
461 /*@C
462    DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems
463 
464    Not Collective
465 
466    Input Argument:
467 .  dm - DM to be used with TS
468 
469    Output Arguments:
470 +  jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
471 -  ctx - context for Jacobian evaluation
472 
473    Level: advanced
474 
475    Note:
476    TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
477    associated with the DM.
478 
479 .seealso: DMTSSetI2Jacobian(),TSGetI2Jacobian()
480 @*/
481 PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
482 {
483   DMTS           tsdm;
484   PetscErrorCode ierr;
485 
486   PetscFunctionBegin;
487   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
488   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
489   if (jac) *jac = tsdm->ops->i2jacobian;
490   if (ctx) *ctx = tsdm->i2jacobianctx;
491   PetscFunctionReturn(0);
492 }
493 
494 /*@C
495    DMTSSetRHSFunction - set TS explicit residual evaluation function
496 
497    Not Collective
498 
499    Input Arguments:
500 +  dm - DM to be used with TS
501 .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
502 -  ctx - context for residual evaluation
503 
504    Level: advanced
505 
506    Note:
507    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
508    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
509    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
510 
511 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
512 @*/
513 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
514 {
515   PetscErrorCode ierr;
516   DMTS           tsdm;
517 
518   PetscFunctionBegin;
519   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
520   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
521   if (func) tsdm->ops->rhsfunction = func;
522   if (ctx)  tsdm->rhsfunctionctx = ctx;
523   PetscFunctionReturn(0);
524 }
525 
526 /*@C
527    DMTSSetRHSFunctionslow - set TS explicit residual evaluation function for slow part
528 
529    Not Collective
530 
531    Input Arguments:
532 +  dm - DM to be used with TS
533 .  func - RHS function evaluation function for slow parts, see TSSetRHSFunctionslow() for calling sequence
534 -  ctx - context for residual evaluation
535 
536    Level: advanced
537 
538    Note:
539    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
540    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
541    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
542 
543 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
544 @*/
545 PetscErrorCode DMTSSetRHSFunctionslow(DM dm,TSRHSFunctionslow func,void *ctx)
546 {
547   PetscErrorCode ierr;
548   DMTS           tsdm;
549 
550   PetscFunctionBegin;
551   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
552   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
553   if (func) tsdm->ops->rhsfunctionslow = func;
554   if (ctx)  tsdm->rhsfunctionslowctx = ctx;
555   PetscFunctionReturn(0);
556 }
557 
558 /*@C
559    DMTSSetRHSFunctionfast - set TS explicit residual evaluation function for fast parts
560 
561    Not Collective
562 
563    Input Arguments:
564 +  dm - DM to be used with TS
565 .  func - RHS function evaluation function for fast parts, see TSSetRHSFunctionfast() for calling sequence
566 -  ctx - context for residual evaluation
567 
568    Level: advanced
569 
570    Note:
571    TSSetRSHFunctionfast() is normally used, but it calls this function internally because the user context is actually
572    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
573    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
574 
575 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
576 @*/
577 PetscErrorCode DMTSSetRHSFunctionfast(DM dm,TSRHSFunctionfast func,void *ctx)
578 {
579   PetscErrorCode ierr;
580   DMTS           tsdm;
581 
582   PetscFunctionBegin;
583   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
584   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
585   if (func) tsdm->ops->rhsfunctionfast = func;
586   if (ctx)  tsdm->rhsfunctionfastctx = ctx;
587   PetscFunctionReturn(0);
588 }
589 
590 /*@C
591    DMTSGetSolutionFunction - gets the TS solution evaluation function
592 
593    Not Collective
594 
595    Input Arguments:
596 .  dm - DM to be used with TS
597 
598    Output Parameters:
599 +  func - solution function evaluation function, see TSSetSolution() for calling sequence
600 -  ctx - context for solution evaluation
601 
602    Level: advanced
603 
604 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSSetSolutionFunction()
605 @*/
606 PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
607 {
608   PetscErrorCode ierr;
609   DMTS           tsdm;
610 
611   PetscFunctionBegin;
612   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
613   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
614   if (func) *func = tsdm->ops->solution;
615   if (ctx)  *ctx  = tsdm->solutionctx;
616   PetscFunctionReturn(0);
617 }
618 
619 /*@C
620    DMTSSetSolutionFunction - set TS solution evaluation function
621 
622    Not Collective
623 
624    Input Arguments:
625 +  dm - DM to be used with TS
626 .  func - solution function evaluation function, see TSSetSolution() for calling sequence
627 -  ctx - context for solution evaluation
628 
629    Level: advanced
630 
631    Note:
632    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
633    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
634    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
635 
636 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSGetSolutionFunction()
637 @*/
638 PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
639 {
640   PetscErrorCode ierr;
641   DMTS           tsdm;
642 
643   PetscFunctionBegin;
644   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
645   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
646   if (func) tsdm->ops->solution = func;
647   if (ctx)  tsdm->solutionctx   = ctx;
648   PetscFunctionReturn(0);
649 }
650 
651 /*@C
652    DMTSSetForcingFunction - set TS forcing function evaluation function
653 
654    Not Collective
655 
656    Input Arguments:
657 +  dm - DM to be used with TS
658 .  f - forcing function evaluation function; see TSForcingFunction
659 -  ctx - context for solution evaluation
660 
661    Level: advanced
662 
663    Note:
664    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
665    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
666    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
667 
668 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
669 @*/
670 PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
671 {
672   PetscErrorCode ierr;
673   DMTS           tsdm;
674 
675   PetscFunctionBegin;
676   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
677   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
678   if (f)    tsdm->ops->forcing = f;
679   if (ctx)  tsdm->forcingctx   = ctx;
680   PetscFunctionReturn(0);
681 }
682 
683 
684 /*@C
685    DMTSGetForcingFunction - get TS forcing function evaluation function
686 
687    Not Collective
688 
689    Input Argument:
690 .   dm - DM to be used with TS
691 
692    Output Arguments:
693 +  f - forcing function evaluation function; see TSForcingFunction for details
694 -  ctx - context for solution evaluation
695 
696    Level: advanced
697 
698    Note:
699    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
700    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
701    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
702 
703 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
704 @*/
705 PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
706 {
707   PetscErrorCode ierr;
708   DMTS           tsdm;
709 
710   PetscFunctionBegin;
711   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
712   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
713   if (f)   *f   = tsdm->ops->forcing;
714   if (ctx) *ctx = tsdm->forcingctx;
715   PetscFunctionReturn(0);
716 }
717 
718 /*@C
719    DMTSGetRHSFunction - get TS explicit residual evaluation function
720 
721    Not Collective
722 
723    Input Argument:
724 .  dm - DM to be used with TS
725 
726    Output Arguments:
727 +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
728 -  ctx - context for residual evaluation
729 
730    Level: advanced
731 
732    Note:
733    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
734    associated with the DM.
735 
736 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
737 @*/
738 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
739 {
740   PetscErrorCode ierr;
741   DMTS           tsdm;
742 
743   PetscFunctionBegin;
744   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
745   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
746   if (func) *func = tsdm->ops->rhsfunction;
747   if (ctx)  *ctx = tsdm->rhsfunctionctx;
748   PetscFunctionReturn(0);
749 }
750 
751 /*@C
752    DMTSGetRHSFunctionslow - get TS explicit residual evaluation function for slow parts
753 
754    Not Collective
755 
756    Input Argument:
757 .  dm - DM to be used with TS
758 
759    Output Arguments:
760 +  func - residual evaluation function, see TSSetRHSFunctionslow() for calling sequence
761 -  ctx - context for residual evaluation
762 
763    Level: advanced
764 
765    Note:
766    TSGetFunctionslow() is normally used, but it calls this function internally because the user context is actually
767    associated with the DM.
768 
769 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
770 @*/
771 PetscErrorCode DMTSGetRHSFunctionslow(DM dm,TSRHSFunctionslow *func,void **ctx)
772 {
773   PetscErrorCode ierr;
774   DMTS           tsdm;
775 
776   PetscFunctionBegin;
777   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
778   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
779   if (func) *func = tsdm->ops->rhsfunctionslow;
780   if (ctx)  *ctx = tsdm->rhsfunctionslowctx;
781   PetscFunctionReturn(0);
782 }
783 
784 /*@C
785    DMTSGetRHSFunctionfast - get TS explicit residual evaluation function for fast part
786 
787    Not Collective
788 
789    Input Argument:
790 .  dm - DM to be used with TS
791 
792    Output Arguments:
793 +  func - residual evaluation function, see TSSetRHSFunctionfast() for calling sequence
794 -  ctx - context for residual evaluation
795 
796    Level: advanced
797 
798    Note:
799    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
800    associated with the DM.
801 
802 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
803 @*/
804 PetscErrorCode DMTSGetRHSFunctionfast(DM dm,TSRHSFunctionfast *func,void **ctx)
805 {
806   PetscErrorCode ierr;
807   DMTS           tsdm;
808 
809   PetscFunctionBegin;
810   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
811   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
812   if (func) *func = tsdm->ops->rhsfunctionfast;
813   if (ctx)  *ctx = tsdm->rhsfunctionfastctx;
814   PetscFunctionReturn(0);
815 }
816 
817 /*@C
818    DMTSSetIJacobian - set TS Jacobian evaluation function
819 
820    Not Collective
821 
822    Input Argument:
823 +  dm - DM to be used with TS
824 .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
825 -  ctx - context for residual evaluation
826 
827    Level: advanced
828 
829    Note:
830    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
831    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
832    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
833 
834 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
835 @*/
836 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
837 {
838   PetscErrorCode ierr;
839   DMTS           sdm;
840 
841   PetscFunctionBegin;
842   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
843   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
844   if (func) sdm->ops->ijacobian = func;
845   if (ctx)  sdm->ijacobianctx   = ctx;
846   PetscFunctionReturn(0);
847 }
848 
849 /*@C
850    DMTSGetIJacobian - get TS Jacobian evaluation function
851 
852    Not Collective
853 
854    Input Argument:
855 .  dm - DM to be used with TS
856 
857    Output Arguments:
858 +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
859 -  ctx - context for residual evaluation
860 
861    Level: advanced
862 
863    Note:
864    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
865    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
866    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
867 
868 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
869 @*/
870 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
871 {
872   PetscErrorCode ierr;
873   DMTS           tsdm;
874 
875   PetscFunctionBegin;
876   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
877   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
878   if (func) *func = tsdm->ops->ijacobian;
879   if (ctx)  *ctx = tsdm->ijacobianctx;
880   PetscFunctionReturn(0);
881 }
882 
883 
884 /*@C
885    DMTSSetRHSJacobian - set TS Jacobian evaluation function
886 
887    Not Collective
888 
889    Input Argument:
890 +  dm - DM to be used with TS
891 .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
892 -  ctx - context for residual evaluation
893 
894    Level: advanced
895 
896    Note:
897    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
898    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
899    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
900 
901 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
902 @*/
903 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
904 {
905   PetscErrorCode ierr;
906   DMTS           tsdm;
907 
908   PetscFunctionBegin;
909   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
910   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
911   if (func) tsdm->ops->rhsjacobian = func;
912   if (ctx)  tsdm->rhsjacobianctx = ctx;
913   PetscFunctionReturn(0);
914 }
915 
916 /*@C
917    DMTSGetRHSJacobian - get TS Jacobian evaluation function
918 
919    Not Collective
920 
921    Input Argument:
922 .  dm - DM to be used with TS
923 
924    Output Arguments:
925 +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
926 -  ctx - context for residual evaluation
927 
928    Level: advanced
929 
930    Note:
931    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
932    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
933    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
934 
935 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
936 @*/
937 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
938 {
939   PetscErrorCode ierr;
940   DMTS           tsdm;
941 
942   PetscFunctionBegin;
943   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
944   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
945   if (func) *func = tsdm->ops->rhsjacobian;
946   if (ctx)  *ctx = tsdm->rhsjacobianctx;
947   PetscFunctionReturn(0);
948 }
949 
950 /*@C
951    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
952 
953    Not Collective
954 
955    Input Arguments:
956 +  dm - DM to be used with TS
957 .  view - viewer function
958 -  load - loading function
959 
960    Level: advanced
961 
962 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
963 @*/
964 PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
965 {
966   PetscErrorCode ierr;
967   DMTS           tsdm;
968 
969   PetscFunctionBegin;
970   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
971   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
972   tsdm->ops->ifunctionview = view;
973   tsdm->ops->ifunctionload = load;
974   PetscFunctionReturn(0);
975 }
976 
977 /*@C
978    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
979 
980    Not Collective
981 
982    Input Arguments:
983 +  dm - DM to be used with TS
984 .  view - viewer function
985 -  load - loading function
986 
987    Level: advanced
988 
989 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
990 @*/
991 PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
992 {
993   PetscErrorCode ierr;
994   DMTS           tsdm;
995 
996   PetscFunctionBegin;
997   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
998   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
999   tsdm->ops->ijacobianview = view;
1000   tsdm->ops->ijacobianload = load;
1001   PetscFunctionReturn(0);
1002 }
1003