1 #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
2 #include <petsc/private/dmimpl.h>
3 
DMTSDestroy(DMTS * kdm)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 = NULL; 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 
DMTSLoad(DMTS kdm,PetscViewer viewer)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 
DMTSView(DMTS kdm,PetscViewer viewer)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);CHKERRQ(ierr);
82     ierr = PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION);CHKERRQ(ierr);
83     ierr = PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION);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);CHKERRQ(ierr);
91     ierr = PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION);CHKERRQ(ierr);
92     ierr = PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION);CHKERRQ(ierr);
93     if (kdm->ops->ijacobianview) {
94       ierr = (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);CHKERRQ(ierr);
95     }
96   }
97   PetscFunctionReturn(0);
98 }
99 
DMTSCreate(MPI_Comm comm,DMTS * kdm)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  */
DMCoarsenHook_DMTS(DM dm,DM dmc,void * ctx)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  */
DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void * ctx)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 
DMSubDomainHook_DMTS(DM dm,DM subdm,void * ctx)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  */
DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void * ctx)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 @*/
DMTSCopy(DMTS kdm,DMTS nkdm)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 @*/
DMGetDMTS(DM dm,DMTS * tsdm)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     (*tsdm)->originaldm = dm;
229     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr);
230     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr);
231   }
232   PetscFunctionReturn(0);
233 }
234 
235 /*@C
236    DMGetDMTSWrite - get write access to private DMTS context from a DM
237 
238    Not Collective
239 
240    Input Argument:
241 .  dm - DM to be used with TS
242 
243    Output Argument:
244 .  tsdm - private DMTS context
245 
246    Level: developer
247 
248 .seealso: DMGetDMTS()
249 @*/
DMGetDMTSWrite(DM dm,DMTS * tsdm)250 PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
251 {
252   PetscErrorCode ierr;
253   DMTS           sdm;
254 
255   PetscFunctionBegin;
256   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
257   ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr);
258   if (!sdm->originaldm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DMTS has a NULL originaldm");
259   if (sdm->originaldm != dm) {  /* Copy on write */
260     DMTS oldsdm = sdm;
261     ierr     = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr);
262     ierr     = DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);CHKERRQ(ierr);
263     ierr     = DMTSCopy(oldsdm,sdm);CHKERRQ(ierr);
264     ierr     = DMTSDestroy((DMTS*)&dm->dmts);CHKERRQ(ierr);
265     dm->dmts = (PetscObject) sdm;
266     sdm->originaldm = dm;
267   }
268   *tsdm = sdm;
269   PetscFunctionReturn(0);
270 }
271 
272 /*@C
273    DMCopyDMTS - copies a DM context to a new DM
274 
275    Logically Collective
276 
277    Input Arguments:
278 +  dmsrc - DM to obtain context from
279 -  dmdest - DM to add context to
280 
281    Level: developer
282 
283    Note:
284    The context is copied by reference. This function does not ensure that a context exists.
285 
286 .seealso: DMGetDMTS(), TSSetDM()
287 @*/
DMCopyDMTS(DM dmsrc,DM dmdest)288 PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
289 {
290   PetscErrorCode ierr;
291 
292   PetscFunctionBegin;
293   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
294   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
295   ierr         = DMTSDestroy((DMTS*)&dmdest->dmts);CHKERRQ(ierr);
296   dmdest->dmts = dmsrc->dmts;
297   ierr         = PetscObjectReference(dmdest->dmts);CHKERRQ(ierr);
298   ierr         = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr);
299   ierr         = DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr);
300   PetscFunctionReturn(0);
301 }
302 
303 /*@C
304    DMTSSetIFunction - set TS implicit function evaluation function
305 
306    Not Collective
307 
308    Input Arguments:
309 +  dm - DM to be used with TS
310 .  func - function evaluating f(t,u,u_t)
311 -  ctx - context for residual evaluation
312 
313    Calling sequence of func:
314 $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
315 
316 +  t   - time at step/stage being solved
317 .  u   - state vector
318 .  u_t - time derivative of state vector
319 .  F   - function vector
320 -  ctx - [optional] user-defined context for matrix evaluation routine
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(), TSSetIFunction(), DMTSSetJacobian()
330 @*/
DMTSSetIFunction(DM dm,TSIFunction func,void * ctx)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 /*@C
345    DMTSGetIFunction - get TS implicit residual evaluation function
346 
347    Not Collective
348 
349    Input Argument:
350 .  dm - DM to be used with TS
351 
352    Output Arguments:
353 +  func - function evaluation function, see TSSetIFunction() for calling sequence
354 -  ctx - context for residual evaluation
355 
356    Level: advanced
357 
358    Note:
359    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
360    associated with the DM.
361 
362 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
363 @*/
DMTSGetIFunction(DM dm,TSIFunction * func,void ** ctx)364 PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
365 {
366   PetscErrorCode ierr;
367   DMTS           tsdm;
368 
369   PetscFunctionBegin;
370   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
371   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
372   if (func) *func = tsdm->ops->ifunction;
373   if (ctx)  *ctx = tsdm->ifunctionctx;
374   PetscFunctionReturn(0);
375 }
376 
377 /*@C
378    DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems
379 
380    Not Collective
381 
382    Input Arguments:
383 +  dm - DM to be used with TS
384 .  fun - function evaluation routine
385 -  ctx - context for residual evaluation
386 
387    Calling sequence of fun:
388 $     PetscErrorCode fun(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,Vec F,ctx);
389 
390 +  t    - time at step/stage being solved
391 .  U    - state vector
392 .  U_t  - time derivative of state vector
393 .  U_tt - second time derivative of state vector
394 .  F    - function vector
395 -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)
396 
397    Level: advanced
398 
399    Note:
400    TSSetI2Function() is normally used, but it calls this function internally because the user context is actually
401    associated with the DM.
402 
403 .seealso: TSSetI2Function()
404 @*/
DMTSSetI2Function(DM dm,TSI2Function fun,void * ctx)405 PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx)
406 {
407   DMTS           tsdm;
408   PetscErrorCode ierr;
409 
410   PetscFunctionBegin;
411   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
412   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
413   if (fun) tsdm->ops->i2function = fun;
414   if (ctx) tsdm->i2functionctx   = ctx;
415   PetscFunctionReturn(0);
416 }
417 
418 /*@C
419    DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems
420 
421    Not Collective
422 
423    Input Argument:
424 .  dm - DM to be used with TS
425 
426    Output Arguments:
427 +  fun - function evaluation function, see TSSetI2Function() for calling sequence
428 -  ctx - context for residual evaluation
429 
430    Level: advanced
431 
432    Note:
433    TSGetI2Function() is normally used, but it calls this function internally because the user context is actually
434    associated with the DM.
435 
436 .seealso: DMTSSetI2Function(),TSGetI2Function()
437 @*/
DMTSGetI2Function(DM dm,TSI2Function * fun,void ** ctx)438 PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx)
439 {
440   DMTS           tsdm;
441   PetscErrorCode ierr;
442 
443   PetscFunctionBegin;
444   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
445   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
446   if (fun) *fun = tsdm->ops->i2function;
447   if (ctx) *ctx = tsdm->i2functionctx;
448   PetscFunctionReturn(0);
449 }
450 
451 /*@C
452    DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems
453 
454    Not Collective
455 
456    Input Arguments:
457 +  dm - DM to be used with TS
458 .  fun - Jacobian evaluation routine
459 -  ctx - context for Jacobian evaluation
460 
461    Calling sequence of jac:
462 $    PetscErrorCode jac(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,PetscReal v,PetscReal a,Mat J,Mat P,void *ctx);
463 
464 +  t    - time at step/stage being solved
465 .  U    - state vector
466 .  U_t  - time derivative of state vector
467 .  U_tt - second time derivative of state vector
468 .  v    - shift for U_t
469 .  a    - shift for U_tt
470 .  J    - Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t  + a*dF/dU_tt
471 .  P    - preconditioning matrix for J, may be same as J
472 -  ctx  - [optional] user-defined context for matrix evaluation routine
473 
474    Level: advanced
475 
476    Note:
477    TSSetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
478    associated with the DM.
479 
480 .seealso: TSSetI2Jacobian()
481 @*/
DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void * ctx)482 PetscErrorCode DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void *ctx)
483 {
484   DMTS           tsdm;
485   PetscErrorCode ierr;
486 
487   PetscFunctionBegin;
488   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
489   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
490   if (jac) tsdm->ops->i2jacobian = jac;
491   if (ctx) tsdm->i2jacobianctx   = ctx;
492   PetscFunctionReturn(0);
493 }
494 
495 /*@C
496    DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems
497 
498    Not Collective
499 
500    Input Argument:
501 .  dm - DM to be used with TS
502 
503    Output Arguments:
504 +  jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
505 -  ctx - context for Jacobian evaluation
506 
507    Level: advanced
508 
509    Note:
510    TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
511    associated with the DM.
512 
513 .seealso: DMTSSetI2Jacobian(),TSGetI2Jacobian()
514 @*/
DMTSGetI2Jacobian(DM dm,TSI2Jacobian * jac,void ** ctx)515 PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
516 {
517   DMTS           tsdm;
518   PetscErrorCode ierr;
519 
520   PetscFunctionBegin;
521   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
522   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
523   if (jac) *jac = tsdm->ops->i2jacobian;
524   if (ctx) *ctx = tsdm->i2jacobianctx;
525   PetscFunctionReturn(0);
526 }
527 
528 /*@C
529    DMTSSetRHSFunction - set TS explicit residual evaluation function
530 
531    Not Collective
532 
533    Input Arguments:
534 +  dm - DM to be used with TS
535 .  func - RHS function evaluation routine
536 -  ctx - context for residual evaluation
537 
538     Calling sequence of func:
539 $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec F,void *ctx);
540 
541 +   ts - timestep context
542 .   t - current timestep
543 .   u - input vector
544 .   F - function vector
545 -   ctx - [optional] user-defined function context
546 
547    Level: advanced
548 
549    Note:
550    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
551    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
552    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
553 
554 .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSSetJacobian()
555 @*/
DMTSSetRHSFunction(DM dm,TSRHSFunction func,void * ctx)556 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
557 {
558   PetscErrorCode ierr;
559   DMTS           tsdm;
560 
561   PetscFunctionBegin;
562   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
563   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
564   if (func) tsdm->ops->rhsfunction = func;
565   if (ctx)  tsdm->rhsfunctionctx = ctx;
566   PetscFunctionReturn(0);
567 }
568 
569 /*@C
570    DMTSSetTransientVariable - sets function to transform from state to transient variables
571 
572    Logically Collective
573 
574    Input Arguments:
575 +  dm - DM to be used with TS
576 .  tvar - a function that transforms to transient variables
577 -  ctx - a context for tvar
578 
579     Calling sequence of tvar:
580 $     PetscErrorCode tvar(TS ts,Vec p,Vec c,void *ctx);
581 
582 +   ts - timestep context
583 .   p - input vector (primative form)
584 .   c - output vector, transient variables (conservative form)
585 -   ctx - [optional] user-defined function context
586 
587    Level: advanced
588 
589    Notes:
590    This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF)
591    can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
592    well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
593    C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
594    evaluated via the chain rule, as in
595 
596      dF/dP + shift * dF/dCdot dC/dP.
597 
598 .seealso: TSSetTransientVariable(), DMTSGetTransientVariable(), DMTSSetIFunction(), DMTSSetIJacobian()
599 @*/
DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void * ctx)600 PetscErrorCode DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void *ctx)
601 {
602   PetscErrorCode ierr;
603   DMTS           dmts;
604 
605   PetscFunctionBegin;
606   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
607   ierr = DMGetDMTSWrite(dm,&dmts);CHKERRQ(ierr);
608   dmts->ops->transientvar = tvar;
609   dmts->transientvarctx = ctx;
610   PetscFunctionReturn(0);
611 }
612 
613 /*@C
614    DMTSGetTransientVariable - gets function to transform from state to transient variables
615 
616    Logically Collective
617 
618    Input Arguments:
619 .  dm - DM to be used with TS
620 
621    Output Arguments:
622 +  tvar - a function that transforms to transient variables
623 -  ctx - a context for tvar
624 
625    Level: advanced
626 
627 .seealso: DMTSSetTransientVariable(), DMTSGetIFunction(), DMTSGetIJacobian()
628 @*/
DMTSGetTransientVariable(DM dm,TSTransientVariable * tvar,void * ctx)629 PetscErrorCode DMTSGetTransientVariable(DM dm,TSTransientVariable *tvar,void *ctx)
630 {
631   PetscErrorCode ierr;
632   DMTS           dmts;
633 
634   PetscFunctionBegin;
635   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
636   ierr = DMGetDMTS(dm,&dmts);CHKERRQ(ierr);
637   if (tvar) *tvar = dmts->ops->transientvar;
638   if (ctx)  *(void**)ctx = dmts->transientvarctx;
639   PetscFunctionReturn(0);
640 }
641 
642 /*@C
643    DMTSGetSolutionFunction - gets the TS solution evaluation function
644 
645    Not Collective
646 
647    Input Arguments:
648 .  dm - DM to be used with TS
649 
650    Output Parameters:
651 +  func - solution function evaluation function, see TSSetSolution() for calling sequence
652 -  ctx - context for solution evaluation
653 
654    Level: advanced
655 
656 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSSetSolutionFunction()
657 @*/
DMTSGetSolutionFunction(DM dm,TSSolutionFunction * func,void ** ctx)658 PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
659 {
660   PetscErrorCode ierr;
661   DMTS           tsdm;
662 
663   PetscFunctionBegin;
664   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
665   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
666   if (func) *func = tsdm->ops->solution;
667   if (ctx)  *ctx  = tsdm->solutionctx;
668   PetscFunctionReturn(0);
669 }
670 
671 /*@C
672    DMTSSetSolutionFunction - set TS solution evaluation function
673 
674    Not Collective
675 
676    Input Arguments:
677 +  dm - DM to be used with TS
678 .  func - solution function evaluation routine
679 -  ctx - context for solution evaluation
680 
681     Calling sequence of f:
682 $     PetscErrorCode f(TS ts,PetscReal t,Vec u,void *ctx);
683 
684 +   ts - timestep context
685 .   t - current timestep
686 .   u - output vector
687 -   ctx - [optional] user-defined function context
688 
689    Level: advanced
690 
691    Note:
692    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
693    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
694    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
695 
696 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSGetSolutionFunction()
697 @*/
DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void * ctx)698 PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
699 {
700   PetscErrorCode ierr;
701   DMTS           tsdm;
702 
703   PetscFunctionBegin;
704   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
705   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
706   if (func) tsdm->ops->solution = func;
707   if (ctx)  tsdm->solutionctx   = ctx;
708   PetscFunctionReturn(0);
709 }
710 
711 /*@C
712    DMTSSetForcingFunction - set TS forcing function evaluation function
713 
714    Not Collective
715 
716    Input Arguments:
717 +  dm - DM to be used with TS
718 .  f - forcing function evaluation routine
719 -  ctx - context for solution evaluation
720 
721     Calling sequence of func:
722 $     PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx);
723 
724 +   ts - timestep context
725 .   t - current timestep
726 .   f - output vector
727 -   ctx - [optional] user-defined function context
728 
729    Level: advanced
730 
731    Note:
732    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
733    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
734    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
735 
736 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
737 @*/
DMTSSetForcingFunction(DM dm,TSForcingFunction f,void * ctx)738 PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
739 {
740   PetscErrorCode ierr;
741   DMTS           tsdm;
742 
743   PetscFunctionBegin;
744   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
745   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
746   if (f)    tsdm->ops->forcing = f;
747   if (ctx)  tsdm->forcingctx   = ctx;
748   PetscFunctionReturn(0);
749 }
750 
751 
752 /*@C
753    DMTSGetForcingFunction - get TS forcing function evaluation function
754 
755    Not Collective
756 
757    Input Argument:
758 .   dm - DM to be used with TS
759 
760    Output Arguments:
761 +  f - forcing function evaluation function; see TSForcingFunction for details
762 -  ctx - context for solution evaluation
763 
764    Level: advanced
765 
766    Note:
767    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
768    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
769    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
770 
771 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
772 @*/
DMTSGetForcingFunction(DM dm,TSForcingFunction * f,void ** ctx)773 PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
774 {
775   PetscErrorCode ierr;
776   DMTS           tsdm;
777 
778   PetscFunctionBegin;
779   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
780   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
781   if (f)   *f   = tsdm->ops->forcing;
782   if (ctx) *ctx = tsdm->forcingctx;
783   PetscFunctionReturn(0);
784 }
785 
786 /*@C
787    DMTSGetRHSFunction - get TS explicit residual evaluation function
788 
789    Not Collective
790 
791    Input Argument:
792 .  dm - DM to be used with TS
793 
794    Output Arguments:
795 +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
796 -  ctx - context for residual evaluation
797 
798    Level: advanced
799 
800    Note:
801    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
802    associated with the DM.
803 
804 .seealso: DMTSSetContext(), DMTSSetRHSFunction(), TSSetRHSFunction()
805 @*/
DMTSGetRHSFunction(DM dm,TSRHSFunction * func,void ** ctx)806 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
807 {
808   PetscErrorCode ierr;
809   DMTS           tsdm;
810 
811   PetscFunctionBegin;
812   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
813   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
814   if (func) *func = tsdm->ops->rhsfunction;
815   if (ctx)  *ctx = tsdm->rhsfunctionctx;
816   PetscFunctionReturn(0);
817 }
818 
819 /*@C
820    DMTSSetIJacobian - set TS Jacobian evaluation function
821 
822    Not Collective
823 
824    Input Argument:
825 +  dm - DM to be used with TS
826 .  func - Jacobian evaluation routine
827 -  ctx - context for residual evaluation
828 
829    Calling sequence of f:
830 $    PetscErrorCode f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx);
831 
832 +  t    - time at step/stage being solved
833 .  U    - state vector
834 .  U_t  - time derivative of state vector
835 .  a    - shift
836 .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
837 .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
838 -  ctx  - [optional] user-defined context for matrix evaluation routine
839 
840    Level: advanced
841 
842    Note:
843    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
844    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
845    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
846 
847 .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSGetJacobian(), TSSetIJacobian(), TSSetIFunction()
848 @*/
DMTSSetIJacobian(DM dm,TSIJacobian func,void * ctx)849 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
850 {
851   PetscErrorCode ierr;
852   DMTS           sdm;
853 
854   PetscFunctionBegin;
855   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
856   ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr);
857   if (func) sdm->ops->ijacobian = func;
858   if (ctx)  sdm->ijacobianctx   = ctx;
859   PetscFunctionReturn(0);
860 }
861 
862 /*@C
863    DMTSGetIJacobian - get TS Jacobian evaluation function
864 
865    Not Collective
866 
867    Input Argument:
868 .  dm - DM to be used with TS
869 
870    Output Arguments:
871 +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
872 -  ctx - context for residual evaluation
873 
874    Level: advanced
875 
876    Note:
877    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
878    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
879    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
880 
881 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
882 @*/
DMTSGetIJacobian(DM dm,TSIJacobian * func,void ** ctx)883 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
884 {
885   PetscErrorCode ierr;
886   DMTS           tsdm;
887 
888   PetscFunctionBegin;
889   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
890   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
891   if (func) *func = tsdm->ops->ijacobian;
892   if (ctx)  *ctx = tsdm->ijacobianctx;
893   PetscFunctionReturn(0);
894 }
895 
896 /*@C
897    DMTSSetRHSJacobian - set TS Jacobian evaluation function
898 
899    Not Collective
900 
901    Input Argument:
902 +  dm - DM to be used with TS
903 .  func - Jacobian evaluation routine
904 -  ctx - context for residual evaluation
905 
906    Calling sequence of func:
907 $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);
908 
909 +  t - current timestep
910 .  u - input vector
911 .  Amat - (approximate) Jacobian matrix
912 .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
913 -  ctx - [optional] user-defined context for matrix evaluation routine
914 
915    Level: advanced
916 
917    Note:
918    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
919    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
920    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
921 
922 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetRHSJacobian()
923 @*/
DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void * ctx)924 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
925 {
926   PetscErrorCode ierr;
927   DMTS           tsdm;
928 
929   PetscFunctionBegin;
930   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
931   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
932   if (func) tsdm->ops->rhsjacobian = func;
933   if (ctx)  tsdm->rhsjacobianctx = ctx;
934   PetscFunctionReturn(0);
935 }
936 
937 /*@C
938    DMTSGetRHSJacobian - get TS Jacobian evaluation function
939 
940    Not Collective
941 
942    Input Argument:
943 .  dm - DM to be used with TS
944 
945    Output Arguments:
946 +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
947 -  ctx - context for residual evaluation
948 
949    Level: advanced
950 
951    Note:
952    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
953    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
954    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
955 
956 .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSSetRHSJacobian(), TSSetRHSJacobian()
957 @*/
DMTSGetRHSJacobian(DM dm,TSRHSJacobian * func,void ** ctx)958 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
959 {
960   PetscErrorCode ierr;
961   DMTS           tsdm;
962 
963   PetscFunctionBegin;
964   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
965   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
966   if (func) *func = tsdm->ops->rhsjacobian;
967   if (ctx)  *ctx = tsdm->rhsjacobianctx;
968   PetscFunctionReturn(0);
969 }
970 
971 /*@C
972    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
973 
974    Not Collective
975 
976    Input Arguments:
977 +  dm - DM to be used with TS
978 .  view - viewer function
979 -  load - loading function
980 
981    Level: advanced
982 
983 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
984 @*/
DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (* view)(void *,PetscViewer),PetscErrorCode (* load)(void **,PetscViewer))985 PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
986 {
987   PetscErrorCode ierr;
988   DMTS           tsdm;
989 
990   PetscFunctionBegin;
991   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
992   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
993   tsdm->ops->ifunctionview = view;
994   tsdm->ops->ifunctionload = load;
995   PetscFunctionReturn(0);
996 }
997 
998 /*@C
999    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
1000 
1001    Not Collective
1002 
1003    Input Arguments:
1004 +  dm - DM to be used with TS
1005 .  view - viewer function
1006 -  load - loading function
1007 
1008    Level: advanced
1009 
1010 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
1011 @*/
DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (* view)(void *,PetscViewer),PetscErrorCode (* load)(void **,PetscViewer))1012 PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
1013 {
1014   PetscErrorCode ierr;
1015   DMTS           tsdm;
1016 
1017   PetscFunctionBegin;
1018   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1019   ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr);
1020   tsdm->ops->ijacobianview = view;
1021   tsdm->ops->ijacobianload = load;
1022   PetscFunctionReturn(0);
1023 }
1024