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