Actual source code: dmts.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1:  #include <petsc/private/tsimpl.h>
  2:  #include <petsc/private/dmimpl.h>

  4: static PetscErrorCode DMTSDestroy(DMTS *kdm)
  5: {

  9:   if (!*kdm) return(0);
 11:   if (--((PetscObject)(*kdm))->refct > 0) {*kdm = 0; return(0);}
 12:   if ((*kdm)->ops->destroy) {((*kdm)->ops->destroy)(*kdm);}
 13:   PetscHeaderDestroy(kdm);
 14:   return(0);
 15: }

 17: PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer)
 18: {

 22:   PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,NULL,PETSC_FUNCTION);
 23:   PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,NULL,PETSC_FUNCTION);
 24:   PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,NULL,PETSC_FUNCTION);
 25:   if (kdm->ops->ifunctionload) {
 26:     (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);
 27:   }
 28:   PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,NULL,PETSC_FUNCTION);
 29:   PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,NULL,PETSC_FUNCTION);
 30:   PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,NULL,PETSC_FUNCTION);
 31:   if (kdm->ops->ijacobianload) {
 32:     (*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer);
 33:   }
 34:   return(0);
 35: }

 37: PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
 38: {
 40:   PetscBool      isascii,isbinary;

 43:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 44:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
 45:   if (isascii) {
 46: #if defined(PETSC_SERIALIZE_FUNCTIONS)
 47:     const char *fname;

 49:     PetscFPTFind(kdm->ops->ifunction,&fname);
 50:     if (fname) {
 51:       PetscViewerASCIIPrintf(viewer,"  IFunction used by TS: %s\n",fname);
 52:     }
 53:     PetscFPTFind(kdm->ops->ijacobian,&fname);
 54:     if (fname) {
 55:       PetscViewerASCIIPrintf(viewer,"  IJacobian function used by TS: %s\n",fname);
 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;

 78:     funcstruct.ifunction         = kdm->ops->ifunction;
 79:     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
 80:     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
 81:     PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION);
 82:     PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION);
 83:     PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION);
 84:     if (kdm->ops->ifunctionview) {
 85:       (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);
 86:     }
 87:     jacstruct.ijacobian = kdm->ops->ijacobian;
 88:     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
 89:     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
 90:     PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION);
 91:     PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION);
 92:     PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION);
 93:     if (kdm->ops->ijacobianview) {
 94:       (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);
 95:     }
 96:   }
 97:   return(0);
 98: }

100: static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
101: {

105:   TSInitializePackage();
106:   PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);
107:   return(0);
108: }

110: /* Attaches the DMTS to the coarse level.
111:  * Under what conditions should we copy versus duplicate?
112:  */
113: static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
114: {

118:   DMCopyDMTS(dm,dmc);
119:   return(0);
120: }

122: /* This could restrict auxiliary information to the coarse level.
123:  */
124: static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
125: {

128:   return(0);
129: }

131: static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
132: {

136:   DMCopyDMTS(dm,subdm);
137:   return(0);
138: }

140: /* This could restrict auxiliary information to the coarse level.
141:  */
142: static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
143: {
145:   return(0);
146: }

148: /*@C
149:    DMTSCopy - copies the information in a DMTS to another DMTS

151:    Not Collective

153:    Input Argument:
154: +  kdm - Original DMTS
155: -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()

157:    Level: developer

159: .seealso: DMTSCreate(), DMTSDestroy()
160: @*/
161: PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
162: {

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;

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;

186:   nkdm->data = kdm->data;

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:   */

194:   /* implementation specific copy hooks */
195:   if (kdm->ops->duplicate) {(*kdm->ops->duplicate)(kdm,nkdm);}
196:   return(0);
197: }

199: /*@C
200:    DMGetDMTS - get read-only private DMTS context from a DM

202:    Not Collective

204:    Input Argument:
205: .  dm - DM to be used with TS

207:    Output Argument:
208: .  tsdm - private DMTS context

210:    Level: developer

212:    Notes:
213:    Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.

215: .seealso: DMGetDMTSWrite()
216: @*/
217: PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
218: {

223:   *tsdm = (DMTS) dm->dmts;
224:   if (!*tsdm) {
225:     PetscInfo(dm,"Creating new DMTS\n");
226:     DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm);
227:     dm->dmts = (PetscObject) *tsdm;
228:     (*tsdm)->originaldm = dm;
229:     DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);
230:     DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);
231:   }
232:   return(0);
233: }

235: /*@C
236:    DMGetDMTSWrite - get write access to private DMTS context from a DM

238:    Not Collective

240:    Input Argument:
241: .  dm - DM to be used with TS

243:    Output Argument:
244: .  tsdm - private DMTS context

246:    Level: developer

248: .seealso: DMGetDMTS()
249: @*/
250: PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
251: {
253:   DMTS           sdm;

257:   DMGetDMTS(dm,&sdm);
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:     PetscInfo(dm,"Copying DMTS due to write\n");
262:     DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);
263:     DMTSCopy(oldsdm,sdm);
264:     DMTSDestroy((DMTS*)&dm->dmts);
265:     dm->dmts = (PetscObject) sdm;
266:     sdm->originaldm = dm;
267:   }
268:   *tsdm = sdm;
269:   return(0);
270: }

272: /*@C
273:    DMCopyDMTS - copies a DM context to a new DM

275:    Logically Collective

277:    Input Arguments:
278: +  dmsrc - DM to obtain context from
279: -  dmdest - DM to add context to

281:    Level: developer

283:    Note:
284:    The context is copied by reference. This function does not ensure that a context exists.

286: .seealso: DMGetDMTS(), TSSetDM()
287: @*/
288: PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
289: {

295:   DMTSDestroy((DMTS*)&dmdest->dmts);
296:   dmdest->dmts = dmsrc->dmts;
297:   PetscObjectReference(dmdest->dmts);
298:   DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);
299:   DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);
300:   return(0);
301: }

303: /*@C
304:    DMTSSetIFunction - set TS implicit function evaluation function

306:    Not Collective

308:    Input Arguments:
309: +  dm - DM to be used with TS
310: .  func - function evaluation function, see TSSetIFunction() for calling sequence
311: -  ctx - context for residual evaluation

313:    Level: advanced

315:    Note:
316:    TSSetFunction() is normally used, but it calls this function internally because the user context is actually
317:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
318:    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.

320: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
321: @*/
322: PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
323: {
325:   DMTS           tsdm;

329:   DMGetDMTSWrite(dm,&tsdm);
330:   if (func) tsdm->ops->ifunction = func;
331:   if (ctx)  tsdm->ifunctionctx = ctx;
332:   return(0);
333: }

335: /*@C
336:    DMTSGetIFunction - get TS implicit residual evaluation function

338:    Not Collective

340:    Input Argument:
341: .  dm - DM to be used with TS

343:    Output Arguments:
344: +  func - function evaluation function, see TSSetIFunction() for calling sequence
345: -  ctx - context for residual evaluation

347:    Level: advanced

349:    Note:
350:    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
351:    associated with the DM.

353: .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
354: @*/
355: PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
356: {
358:   DMTS           tsdm;

362:   DMGetDMTS(dm,&tsdm);
363:   if (func) *func = tsdm->ops->ifunction;
364:   if (ctx)  *ctx = tsdm->ifunctionctx;
365:   return(0);
366: }

368: /*@C
369:    DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems

371:    Not Collective

373:    Input Arguments:
374: +  dm - DM to be used with TS
375: .  fun - function evaluation function, see TSSetI2Function() for calling sequence
376: -  ctx - context for residual evaluation

378:    Level: advanced

380:    Note:
381:    TSSetI2Function() is normally used, but it calls this function internally because the user context is actually
382:    associated with the DM.

384: .seealso: TSSetI2Function()
385: @*/
386: PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx)
387: {
388:   DMTS           tsdm;

393:   DMGetDMTSWrite(dm,&tsdm);
394:   if (fun) tsdm->ops->i2function = fun;
395:   if (ctx) tsdm->i2functionctx   = ctx;
396:   return(0);
397: }

399: /*@C
400:    DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems

402:    Not Collective

404:    Input Argument:
405: .  dm - DM to be used with TS

407:    Output Arguments:
408: +  fun - function evaluation function, see TSSetI2Function() for calling sequence
409: -  ctx - context for residual evaluation

411:    Level: advanced

413:    Note:
414:    TSGetI2Function() is normally used, but it calls this function internally because the user context is actually
415:    associated with the DM.

417: .seealso: DMTSSetI2Function(),TSGetI2Function()
418: @*/
419: PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx)
420: {
421:   DMTS           tsdm;

426:   DMGetDMTS(dm,&tsdm);
427:   if (fun) *fun = tsdm->ops->i2function;
428:   if (ctx) *ctx = tsdm->i2functionctx;
429:   return(0);
430: }

432: /*@C
433:    DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems

435:    Not Collective

437:    Input Arguments:
438: +  dm - DM to be used with TS
439: .  fun - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
440: -  ctx - context for Jacobian evaluation

442:    Level: advanced

444:    Note:
445:    TSSetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
446:    associated with the DM.

448: .seealso: TSSetI2Jacobian()
449: @*/
450: PetscErrorCode DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void *ctx)
451: {
452:   DMTS           tsdm;

457:   DMGetDMTSWrite(dm,&tsdm);
458:   if (jac) tsdm->ops->i2jacobian = jac;
459:   if (ctx) tsdm->i2jacobianctx   = ctx;
460:   return(0);
461: }

463: /*@C
464:    DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems

466:    Not Collective

468:    Input Argument:
469: .  dm - DM to be used with TS

471:    Output Arguments:
472: +  jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
473: -  ctx - context for Jacobian evaluation

475:    Level: advanced

477:    Note:
478:    TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
479:    associated with the DM.

481: .seealso: DMTSSetI2Jacobian(),TSGetI2Jacobian()
482: @*/
483: PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
484: {
485:   DMTS           tsdm;

490:   DMGetDMTS(dm,&tsdm);
491:   if (jac) *jac = tsdm->ops->i2jacobian;
492:   if (ctx) *ctx = tsdm->i2jacobianctx;
493:   return(0);
494: }

496: /*@C
497:    DMTSSetRHSFunction - set TS explicit residual evaluation function

499:    Not Collective

501:    Input Arguments:
502: +  dm - DM to be used with TS
503: .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
504: -  ctx - context for residual evaluation

506:    Level: advanced

508:    Note:
509:    TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
510:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
511:    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.

513: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
514: @*/
515: PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
516: {
518:   DMTS           tsdm;

522:   DMGetDMTSWrite(dm,&tsdm);
523:   if (func) tsdm->ops->rhsfunction = func;
524:   if (ctx)  tsdm->rhsfunctionctx = ctx;
525:   return(0);
526: }

528: /*@C
529:    DMTSSetTransientVariable - sets function to transform from state to transient variables

531:    Logically Collective

533:    Input Arguments:
534: +  dm - DM to be used with TS
535: .  tvar - a function that transforms in-place to transient variables
536: -  ctx - a context for tvar

538:    Level: advanced

540:    Notes:
541:    This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF)
542:    can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
543:    well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
544:    C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
545:    evaluated via the chain rule, as in

547:      dF/dP + shift * dF/dCdot dC/dP.

549: .seealso: TSSetTransientVariable(), DMTSGetTransientVariable(), DMTSSetIFunction(), DMTSSetIJacobian()
550: @*/
551: PetscErrorCode DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void *ctx)
552: {
554:   DMTS           dmts;

558:   DMGetDMTSWrite(dm,&dmts);
559:   dmts->ops->transientvar = tvar;
560:   dmts->transientvarctx = ctx;
561:   return(0);
562: }

564: /*@C
565:    DMTSGetTransientVariable - gets function to transform from state to transient variables

567:    Logically Collective

569:    Input Arguments:
570: .  dm - DM to be used with TS

572:    Output Arguments:
573: +  tvar - a function that transforms in-place to transient variables
574: -  ctx - a context for tvar

576:    Level: advanced

578: .seealso: DMTSSetTransientVariable(), DMTSGetIFunction(), DMTSGetIJacobian()
579: @*/
580: PetscErrorCode DMTSGetTransientVariable(DM dm,TSTransientVariable *tvar,void *ctx)
581: {
583:   DMTS           dmts;

587:   DMGetDMTS(dm,&dmts);
588:   if (tvar) *tvar = dmts->ops->transientvar;
589:   if (ctx)  *(void**)ctx = dmts->transientvarctx;
590:   return(0);
591: }

593: /*@C
594:    DMTSGetSolutionFunction - gets the TS solution evaluation function

596:    Not Collective

598:    Input Arguments:
599: .  dm - DM to be used with TS

601:    Output Parameters:
602: +  func - solution function evaluation function, see TSSetSolution() for calling sequence
603: -  ctx - context for solution evaluation

605:    Level: advanced

607: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSSetSolutionFunction()
608: @*/
609: PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
610: {
612:   DMTS           tsdm;

616:   DMGetDMTS(dm,&tsdm);
617:   if (func) *func = tsdm->ops->solution;
618:   if (ctx)  *ctx  = tsdm->solutionctx;
619:   return(0);
620: }

622: /*@C
623:    DMTSSetSolutionFunction - set TS solution evaluation function

625:    Not Collective

627:    Input Arguments:
628: +  dm - DM to be used with TS
629: .  func - solution function evaluation function, see TSSetSolution() for calling sequence
630: -  ctx - context for solution evaluation

632:    Level: advanced

634:    Note:
635:    TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
636:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
637:    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.

639: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSGetSolutionFunction()
640: @*/
641: PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
642: {
644:   DMTS           tsdm;

648:   DMGetDMTSWrite(dm,&tsdm);
649:   if (func) tsdm->ops->solution = func;
650:   if (ctx)  tsdm->solutionctx   = ctx;
651:   return(0);
652: }

654: /*@C
655:    DMTSSetForcingFunction - set TS forcing function evaluation function

657:    Not Collective

659:    Input Arguments:
660: +  dm - DM to be used with TS
661: .  f - forcing function evaluation function; see TSForcingFunction
662: -  ctx - context for solution evaluation

664:    Level: advanced

666:    Note:
667:    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
668:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
669:    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.

671: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
672: @*/
673: PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
674: {
676:   DMTS           tsdm;

680:   DMGetDMTSWrite(dm,&tsdm);
681:   if (f)    tsdm->ops->forcing = f;
682:   if (ctx)  tsdm->forcingctx   = ctx;
683:   return(0);
684: }


687: /*@C
688:    DMTSGetForcingFunction - get TS forcing function evaluation function

690:    Not Collective

692:    Input Argument:
693: .   dm - DM to be used with TS

695:    Output Arguments:
696: +  f - forcing function evaluation function; see TSForcingFunction for details
697: -  ctx - context for solution evaluation

699:    Level: advanced

701:    Note:
702:    TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
703:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
704:    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.

706: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
707: @*/
708: PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
709: {
711:   DMTS           tsdm;

715:   DMGetDMTSWrite(dm,&tsdm);
716:   if (f)   *f   = tsdm->ops->forcing;
717:   if (ctx) *ctx = tsdm->forcingctx;
718:   return(0);
719: }

721: /*@C
722:    DMTSGetRHSFunction - get TS explicit residual evaluation function

724:    Not Collective

726:    Input Argument:
727: .  dm - DM to be used with TS

729:    Output Arguments:
730: +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
731: -  ctx - context for residual evaluation

733:    Level: advanced

735:    Note:
736:    TSGetFunction() is normally used, but it calls this function internally because the user context is actually
737:    associated with the DM.

739: .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
740: @*/
741: PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
742: {
744:   DMTS           tsdm;

748:   DMGetDMTS(dm,&tsdm);
749:   if (func) *func = tsdm->ops->rhsfunction;
750:   if (ctx)  *ctx = tsdm->rhsfunctionctx;
751:   return(0);
752: }

754: /*@C
755:    DMTSSetIJacobian - set TS Jacobian evaluation function

757:    Not Collective

759:    Input Argument:
760: +  dm - DM to be used with TS
761: .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
762: -  ctx - context for residual evaluation

764:    Level: advanced

766:    Note:
767:    TSSetJacobian() 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 Jacobian.

771: .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
772: @*/
773: PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
774: {
776:   DMTS           sdm;

780:   DMGetDMTSWrite(dm,&sdm);
781:   if (func) sdm->ops->ijacobian = func;
782:   if (ctx)  sdm->ijacobianctx   = ctx;
783:   return(0);
784: }

786: /*@C
787:    DMTSGetIJacobian - get TS Jacobian evaluation function

789:    Not Collective

791:    Input Argument:
792: .  dm - DM to be used with TS

794:    Output Arguments:
795: +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
796: -  ctx - context for residual evaluation

798:    Level: advanced

800:    Note:
801:    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
802:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
803:    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.

805: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
806: @*/
807: PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
808: {
810:   DMTS           tsdm;

814:   DMGetDMTS(dm,&tsdm);
815:   if (func) *func = tsdm->ops->ijacobian;
816:   if (ctx)  *ctx = tsdm->ijacobianctx;
817:   return(0);
818: }


821: /*@C
822:    DMTSSetRHSJacobian - set TS Jacobian evaluation function

824:    Not Collective

826:    Input Argument:
827: +  dm - DM to be used with TS
828: .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
829: -  ctx - context for residual evaluation

831:    Level: advanced

833:    Note:
834:    TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
835:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
836:    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.

838: .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
839: @*/
840: PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
841: {
843:   DMTS           tsdm;

847:   DMGetDMTSWrite(dm,&tsdm);
848:   if (func) tsdm->ops->rhsjacobian = func;
849:   if (ctx)  tsdm->rhsjacobianctx = ctx;
850:   return(0);
851: }

853: /*@C
854:    DMTSGetRHSJacobian - get TS Jacobian evaluation function

856:    Not Collective

858:    Input Argument:
859: .  dm - DM to be used with TS

861:    Output Arguments:
862: +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
863: -  ctx - context for residual evaluation

865:    Level: advanced

867:    Note:
868:    TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
869:    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
870:    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.

872: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
873: @*/
874: PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
875: {
877:   DMTS           tsdm;

881:   DMGetDMTS(dm,&tsdm);
882:   if (func) *func = tsdm->ops->rhsjacobian;
883:   if (ctx)  *ctx = tsdm->rhsjacobianctx;
884:   return(0);
885: }

887: /*@C
888:    DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context

890:    Not Collective

892:    Input Arguments:
893: +  dm - DM to be used with TS
894: .  view - viewer function
895: -  load - loading function

897:    Level: advanced

899: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
900: @*/
901: PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
902: {
904:   DMTS           tsdm;

908:   DMGetDMTSWrite(dm,&tsdm);
909:   tsdm->ops->ifunctionview = view;
910:   tsdm->ops->ifunctionload = load;
911:   return(0);
912: }

914: /*@C
915:    DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context

917:    Not Collective

919:    Input Arguments:
920: +  dm - DM to be used with TS
921: .  view - viewer function
922: -  load - loading function

924:    Level: advanced

926: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
927: @*/
928: PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
929: {
931:   DMTS           tsdm;

935:   DMGetDMTSWrite(dm,&tsdm);
936:   tsdm->ops->ijacobianview = view;
937:   tsdm->ops->ijacobianload = load;
938:   return(0);
939: }