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