xref: /petsc/src/ts/utils/dmts.c (revision f8e4bde84cdee892096ab78baf4256ae384d8e42)
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   PetscInt       num = 1;
24   PetscErrorCode ierr;
25 
26   PetscFunctionBegin;
27   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,&num,PETSC_FUNCTION);CHKERRQ(ierr);
28   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,&num,PETSC_FUNCTION);CHKERRQ(ierr);
29   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,&num,PETSC_FUNCTION);CHKERRQ(ierr);
30   if (kdm->ops->ifunctionload) {
31     ierr = (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);CHKERRQ(ierr);
32   }
33   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,&num,PETSC_FUNCTION);CHKERRQ(ierr);
34   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,&num,PETSC_FUNCTION);CHKERRQ(ierr);
35   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,&num,PETSC_FUNCTION);CHKERRQ(ierr);
36   if (kdm->ops->ijacobianload) {
37     ierr = (*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer);CHKERRQ(ierr);
38   }
39   PetscFunctionReturn(0);
40 }
41 
42 #undef __FUNCT__
43 #define __FUNCT__ "DMTSView"
44 PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
45 {
46   PetscErrorCode ierr;
47   PetscBool      isascii,isbinary;
48 
49   PetscFunctionBegin;
50   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
51   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
52   if (isascii) {
53 #if defined(PETSC_SERIALIZE_FUNCTIONS)
54     const char *fname;
55 
56     ierr = PetscFPTFind(kdm->ops->ifunction,&fname);CHKERRQ(ierr);
57     if (fname) {
58       ierr = PetscViewerASCIIPrintf(viewer,"  IFunction used by TS: %s\n",fname);CHKERRQ(ierr);
59     }
60     ierr = PetscFPTFind(kdm->ops->ijacobian,&fname);CHKERRQ(ierr);
61     if (fname) {
62       ierr = PetscViewerASCIIPrintf(viewer,"  IJacobian function used by TS: %s\n",fname);CHKERRQ(ierr);
63     }
64 #endif
65   } else if (isbinary) {
66     struct {
67       TSIFunction ifunction;
68     } funcstruct;
69     struct {
70       PetscErrorCode (*ifunctionview)(void*,PetscViewer);
71     } funcviewstruct;
72     struct {
73       PetscErrorCode (*ifunctionload)(void**,PetscViewer);
74     } funcloadstruct;
75     struct {
76       TSIJacobian ijacobian;
77     } jacstruct;
78     struct {
79       PetscErrorCode (*ijacobianview)(void*,PetscViewer);
80     } jacviewstruct;
81     struct {
82       PetscErrorCode (*ijacobianload)(void**,PetscViewer);
83     } jacloadstruct;
84 
85     funcstruct.ifunction         = kdm->ops->ifunction;
86     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
87     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
88     ierr = PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
89     ierr = PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
90     ierr = PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
91     if (kdm->ops->ifunctionview) {
92       ierr = (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);CHKERRQ(ierr);
93     }
94     jacstruct.ijacobian = kdm->ops->ijacobian;
95     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
96     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
97     ierr = PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
98     ierr = PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
99     ierr = PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
100     if (kdm->ops->ijacobianview) {
101       ierr = (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);CHKERRQ(ierr);
102     }
103   }
104   PetscFunctionReturn(0);
105 }
106 
107 #undef __FUNCT__
108 #define __FUNCT__ "DMTSCreate"
109 static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
110 {
111   PetscErrorCode ierr;
112 
113   PetscFunctionBegin;
114   ierr = TSInitializePackage();CHKERRQ(ierr);
115   ierr = PetscHeaderCreate(*kdm, _p_DMTS, struct _DMTSOps, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);CHKERRQ(ierr);
116   ierr = PetscMemzero((*kdm)->ops, sizeof(struct _DMTSOps));CHKERRQ(ierr);
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__ "DMTSSetIFunction"
325 /*@C
326    DMTSSetIFunction - set TS implicit function evaluation function
327 
328    Not Collective
329 
330    Input Arguments:
331 +  dm - DM to be used with TS
332 .  func - function evaluation function, see TSSetIFunction() for calling sequence
333 -  ctx - context for residual evaluation
334 
335    Level: advanced
336 
337    Note:
338    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
339    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
340    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
341 
342 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
343 @*/
344 PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
345 {
346   PetscErrorCode ierr;
347   DMTS           tsdm;
348 
349   PetscFunctionBegin;
350   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
351   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
352   if (func) tsdm->ops->ifunction = func;
353   if (ctx)  tsdm->ifunctionctx = ctx;
354   PetscFunctionReturn(0);
355 }
356 
357 #undef __FUNCT__
358 #define __FUNCT__ "DMTSGetIFunction"
359 /*@C
360    DMTSGetIFunction - get TS implicit residual evaluation function
361 
362    Not Collective
363 
364    Input Argument:
365 .  dm - DM to be used with TS
366 
367    Output Arguments:
368 +  func - function evaluation function, see TSSetIFunction() for calling sequence
369 -  ctx - context for residual evaluation
370 
371    Level: advanced
372 
373    Note:
374    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
375    associated with the DM.
376 
377 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
378 @*/
379 PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
380 {
381   PetscErrorCode ierr;
382   DMTS           tsdm;
383 
384   PetscFunctionBegin;
385   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
386   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
387   if (func) *func = tsdm->ops->ifunction;
388   if (ctx)  *ctx = tsdm->ifunctionctx;
389   PetscFunctionReturn(0);
390 }
391 
392 
393 #undef __FUNCT__
394 #define __FUNCT__ "DMTSSetRHSFunction"
395 /*@C
396    DMTSSetRHSFunction - set TS explicit residual evaluation function
397 
398    Not Collective
399 
400    Input Arguments:
401 +  dm - DM to be used with TS
402 .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
403 -  ctx - context for residual evaluation
404 
405    Level: advanced
406 
407    Note:
408    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
409    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
410    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
411 
412 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
413 @*/
414 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
415 {
416   PetscErrorCode ierr;
417   DMTS           tsdm;
418 
419   PetscFunctionBegin;
420   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
421   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
422   if (func) tsdm->ops->rhsfunction = func;
423   if (ctx)  tsdm->rhsfunctionctx = ctx;
424   PetscFunctionReturn(0);
425 }
426 
427 #undef __FUNCT__
428 #define __FUNCT__ "DMTSGetSolutionFunction"
429 /*@C
430    DMTSGetSolutionFunction - gets the TS solution evaluation function
431 
432    Not Collective
433 
434    Input Arguments:
435 .  dm - DM to be used with TS
436 
437    Output Parameters:
438 +  func - solution function evaluation function, see TSSetSolution() for calling sequence
439 -  ctx - context for solution evaluation
440 
441    Level: advanced
442 
443 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
444 @*/
445 PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
446 {
447   PetscErrorCode ierr;
448   DMTS           tsdm;
449 
450   PetscFunctionBegin;
451   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
452   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
453   if (func) *func = tsdm->ops->solution;
454   if (ctx)  *ctx  = tsdm->solutionctx;
455   PetscFunctionReturn(0);
456 }
457 
458 #undef __FUNCT__
459 #define __FUNCT__ "DMTSSetSolutionFunction"
460 /*@C
461    DMTSSetSolutionFunction - set TS solution evaluation function
462 
463    Not Collective
464 
465    Input Arguments:
466 +  dm - DM to be used with TS
467 .  func - solution function evaluation function, see TSSetSolution() for calling sequence
468 -  ctx - context for solution evaluation
469 
470    Level: advanced
471 
472    Note:
473    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
474    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
475    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
476 
477 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
478 @*/
479 PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
480 {
481   PetscErrorCode ierr;
482   DMTS           tsdm;
483 
484   PetscFunctionBegin;
485   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
486   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
487   if (func) tsdm->ops->solution = func;
488   if (ctx)  tsdm->solutionctx   = ctx;
489   PetscFunctionReturn(0);
490 }
491 
492 #undef __FUNCT__
493 #define __FUNCT__ "DMTSSetForcingFunction"
494 /*@C
495    DMTSSetForcingFunction - set TS forcing function evaluation function
496 
497    Not Collective
498 
499    Input Arguments:
500 +  dm - DM to be used with TS
501 .  f - forcing function evaluation function; see TSForcingFunction
502 -  ctx - context for solution evaluation
503 
504    Level: advanced
505 
506    Note:
507    TSSetForcingFunction() 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(), TSSetForcingFunction(), DMTSGetForcingFunction()
512 @*/
513 PetscErrorCode DMTSSetForcingFunction(DM dm,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),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 (f) tsdm->ops->forcing = f;
522   if (ctx)  tsdm->forcingctx   = ctx;
523   PetscFunctionReturn(0);
524 }
525 
526 
527 #undef __FUNCT__
528 #define __FUNCT__ "DMTSGetForcingFunction"
529 /*@C
530    DMTSGetForcingFunction - get TS forcing function evaluation function
531 
532    Not Collective
533 
534    Input Argument:
535 .   dm - DM to be used with TS
536 
537    Output Arguments:
538 +  f - forcing function evaluation function; see TSForcingFunction for details
539 -  ctx - context for solution evaluation
540 
541    Level: advanced
542 
543    Note:
544    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
545    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
546    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
547 
548 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
549 @*/
550 PetscErrorCode DMTSGetForcingFunction(DM dm,PetscErrorCode (**f)(TS,PetscReal,Vec,void*),void **ctx)
551 {
552   PetscErrorCode ierr;
553   DMTS           tsdm;
554 
555   PetscFunctionBegin;
556   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
557   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
558   if (f) *f = tsdm->ops->forcing;
559   if (ctx) *ctx = tsdm->forcingctx;
560   PetscFunctionReturn(0);
561 }
562 
563 #undef __FUNCT__
564 #define __FUNCT__ "DMTSGetRHSFunction"
565 /*@C
566    DMTSGetRHSFunction - get TS explicit residual evaluation function
567 
568    Not Collective
569 
570    Input Argument:
571 .  dm - DM to be used with TS
572 
573    Output Arguments:
574 +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
575 -  ctx - context for residual evaluation
576 
577    Level: advanced
578 
579    Note:
580    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
581    associated with the DM.
582 
583 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
584 @*/
585 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
586 {
587   PetscErrorCode ierr;
588   DMTS           tsdm;
589 
590   PetscFunctionBegin;
591   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
592   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
593   if (func) *func = tsdm->ops->rhsfunction;
594   if (ctx)  *ctx = tsdm->rhsfunctionctx;
595   PetscFunctionReturn(0);
596 }
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "DMTSSetIJacobian"
600 /*@C
601    DMTSSetIJacobian - set TS Jacobian evaluation function
602 
603    Not Collective
604 
605    Input Argument:
606 +  dm - DM to be used with TS
607 .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
608 -  ctx - context for residual evaluation
609 
610    Level: advanced
611 
612    Note:
613    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
614    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
615    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
616 
617 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
618 @*/
619 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
620 {
621   PetscErrorCode ierr;
622   DMTS           sdm;
623 
624   PetscFunctionBegin;
625   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
626   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
627   if (func) sdm->ops->ijacobian = func;
628   if (ctx)  sdm->ijacobianctx   = ctx;
629   PetscFunctionReturn(0);
630 }
631 
632 #undef __FUNCT__
633 #define __FUNCT__ "DMTSGetIJacobian"
634 /*@C
635    DMTSGetIJacobian - get TS Jacobian evaluation function
636 
637    Not Collective
638 
639    Input Argument:
640 .  dm - DM to be used with TS
641 
642    Output Arguments:
643 +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
644 -  ctx - context for residual evaluation
645 
646    Level: advanced
647 
648    Note:
649    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
650    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
651    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
652 
653 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
654 @*/
655 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
656 {
657   PetscErrorCode ierr;
658   DMTS           tsdm;
659 
660   PetscFunctionBegin;
661   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
662   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
663   if (func) *func = tsdm->ops->ijacobian;
664   if (ctx)  *ctx = tsdm->ijacobianctx;
665   PetscFunctionReturn(0);
666 }
667 
668 
669 #undef __FUNCT__
670 #define __FUNCT__ "DMTSSetRHSJacobian"
671 /*@C
672    DMTSSetRHSJacobian - set TS Jacobian evaluation function
673 
674    Not Collective
675 
676    Input Argument:
677 +  dm - DM to be used with TS
678 .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
679 -  ctx - context for residual evaluation
680 
681    Level: advanced
682 
683    Note:
684    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
685    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
686    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
687 
688 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
689 @*/
690 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
691 {
692   PetscErrorCode ierr;
693   DMTS           tsdm;
694 
695   PetscFunctionBegin;
696   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
697   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
698   if (func) tsdm->ops->rhsjacobian = func;
699   if (ctx)  tsdm->rhsjacobianctx = ctx;
700   PetscFunctionReturn(0);
701 }
702 
703 #undef __FUNCT__
704 #define __FUNCT__ "DMTSGetRHSJacobian"
705 /*@C
706    DMTSGetRHSJacobian - get TS Jacobian evaluation function
707 
708    Not Collective
709 
710    Input Argument:
711 .  dm - DM to be used with TS
712 
713    Output Arguments:
714 +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
715 -  ctx - context for residual evaluation
716 
717    Level: advanced
718 
719    Note:
720    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
721    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
722    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
723 
724 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
725 @*/
726 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
727 {
728   PetscErrorCode ierr;
729   DMTS           tsdm;
730 
731   PetscFunctionBegin;
732   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
733   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
734   if (func) *func = tsdm->ops->rhsjacobian;
735   if (ctx)  *ctx = tsdm->rhsjacobianctx;
736   PetscFunctionReturn(0);
737 }
738 
739 #undef __FUNCT__
740 #define __FUNCT__ "DMTSSetIFunctionSerialize"
741 /*@C
742    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
743 
744    Not Collective
745 
746    Input Arguments:
747 +  dm - DM to be used with TS
748 .  view - viewer function
749 -  load - loading function
750 
751    Level: advanced
752 
753 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
754 @*/
755 PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
756 {
757   PetscErrorCode ierr;
758   DMTS           tsdm;
759 
760   PetscFunctionBegin;
761   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
762   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
763   tsdm->ops->ifunctionview = view;
764   tsdm->ops->ifunctionload = load;
765   PetscFunctionReturn(0);
766 }
767 
768 #undef __FUNCT__
769 #define __FUNCT__ "DMTSSetIJacobianSerialize"
770 /*@C
771    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
772 
773    Not Collective
774 
775    Input Arguments:
776 +  dm - DM to be used with TS
777 .  view - viewer function
778 -  load - loading function
779 
780    Level: advanced
781 
782 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
783 @*/
784 PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
785 {
786   PetscErrorCode ierr;
787   DMTS           tsdm;
788 
789   PetscFunctionBegin;
790   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
791   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
792   tsdm->ops->ijacobianview = view;
793   tsdm->ops->ijacobianload = load;
794   PetscFunctionReturn(0);
795 }
796