Actual source code: dmts.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  1: #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
  2: #include <petsc/private/dmimpl.h>

  6: static PetscErrorCode DMTSDestroy(DMTS *kdm)
  7: {

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

 21: PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer)
 22: {

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

 43: PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
 44: {
 46:   PetscBool      isascii,isbinary;

 49:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 50:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
 51:   if (isascii) {
 52: #if defined(PETSC_SERIALIZE_FUNCTIONS)
 53:     const char *fname;

 55:     PetscFPTFind(kdm->ops->ifunction,&fname);
 56:     if (fname) {
 57:       PetscViewerASCIIPrintf(viewer,"  IFunction used by TS: %s\n",fname);
 58:     }
 59:     PetscFPTFind(kdm->ops->ijacobian,&fname);
 60:     if (fname) {
 61:       PetscViewerASCIIPrintf(viewer,"  IJacobian function used by TS: %s\n",fname);
 62:     }
 63: #endif
 64:   } else if (isbinary) {
 65:     struct {
 66:       TSIFunction ifunction;
 67:     } funcstruct;
 68:     struct {
 69:       PetscErrorCode (*ifunctionview)(void*,PetscViewer);
 70:     } funcviewstruct;
 71:     struct {
 72:       PetscErrorCode (*ifunctionload)(void**,PetscViewer);
 73:     } funcloadstruct;
 74:     struct {
 75:       TSIJacobian ijacobian;
 76:     } jacstruct;
 77:     struct {
 78:       PetscErrorCode (*ijacobianview)(void*,PetscViewer);
 79:     } jacviewstruct;
 80:     struct {
 81:       PetscErrorCode (*ijacobianload)(void**,PetscViewer);
 82:     } jacloadstruct;

 84:     funcstruct.ifunction         = kdm->ops->ifunction;
 85:     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
 86:     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
 87:     PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION,PETSC_FALSE);
 88:     PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION,PETSC_FALSE);
 89:     PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION,PETSC_FALSE);
 90:     if (kdm->ops->ifunctionview) {
 91:       (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);
 92:     }
 93:     jacstruct.ijacobian = kdm->ops->ijacobian;
 94:     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
 95:     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
 96:     PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION,PETSC_FALSE);
 97:     PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION,PETSC_FALSE);
 98:     PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION,PETSC_FALSE);
 99:     if (kdm->ops->ijacobianview) {
100:       (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);
101:     }
102:   }
103:   return(0);
104: }

108: static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
109: {

113:   TSInitializePackage();
114:   PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);
115:   return(0);
116: }

120: /* Attaches the DMTS to the coarse level.
121:  * Under what conditions should we copy versus duplicate?
122:  */
123: static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
124: {

128:   DMCopyDMTS(dm,dmc);
129:   return(0);
130: }

134: /* This could restrict auxiliary information to the coarse level.
135:  */
136: static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
137: {

140:   return(0);
141: }

145: static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
146: {

150:   DMCopyDMTS(dm,subdm);
151:   return(0);
152: }

156: /* This could restrict auxiliary information to the coarse level.
157:  */
158: static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
159: {
161:   return(0);
162: }

166: /*@C
167:    DMTSCopy - copies the information in a DMTS to another DMTS

169:    Not Collective

171:    Input Argument:
172: +  kdm - Original DMTS
173: -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()

175:    Level: developer

177: .seealso: DMTSCreate(), DMTSDestroy()
178: @*/
179: PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
180: {

186:   nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
187:   nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
188:   nkdm->ops->ifunction   = kdm->ops->ifunction;
189:   nkdm->ops->ijacobian   = kdm->ops->ijacobian;
190:   nkdm->ops->i2function  = kdm->ops->i2function;
191:   nkdm->ops->i2jacobian  = kdm->ops->i2jacobian;
192:   nkdm->ops->solution    = kdm->ops->solution;
193:   nkdm->ops->destroy     = kdm->ops->destroy;
194:   nkdm->ops->duplicate   = kdm->ops->duplicate;

196:   nkdm->rhsfunctionctx = kdm->rhsfunctionctx;
197:   nkdm->rhsjacobianctx = kdm->rhsjacobianctx;
198:   nkdm->ifunctionctx   = kdm->ifunctionctx;
199:   nkdm->ijacobianctx   = kdm->ijacobianctx;
200:   nkdm->i2functionctx  = kdm->i2functionctx;
201:   nkdm->i2jacobianctx  = kdm->i2jacobianctx;
202:   nkdm->solutionctx    = kdm->solutionctx;

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

206:   /*
207:   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
208:   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
209:   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
210:   */

212:   /* implementation specific copy hooks */
213:   if (kdm->ops->duplicate) {(*kdm->ops->duplicate)(kdm,nkdm);}
214:   return(0);
215: }

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

222:    Not Collective

224:    Input Argument:
225: .  dm - DM to be used with TS

227:    Output Argument:
228: .  tsdm - private DMTS context

230:    Level: developer

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

235: .seealso: DMGetDMTSWrite()
236: @*/
237: PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
238: {

243:   *tsdm = (DMTS) dm->dmts;
244:   if (!*tsdm) {
245:     PetscInfo(dm,"Creating new DMTS\n");
246:     DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm);
247:     dm->dmts = (PetscObject) *tsdm;
248:     DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);
249:     DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);
250:   }
251:   return(0);
252: }

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

259:    Not Collective

261:    Input Argument:
262: .  dm - DM to be used with TS

264:    Output Argument:
265: .  tsdm - private DMTS context

267:    Level: developer

269: .seealso: DMGetDMTS()
270: @*/
271: PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
272: {
274:   DMTS           sdm;

278:   DMGetDMTS(dm,&sdm);
279:   if (!sdm->originaldm) sdm->originaldm = dm;
280:   if (sdm->originaldm != dm) {  /* Copy on write */
281:     DMTS oldsdm = sdm;
282:     PetscInfo(dm,"Copying DMTS due to write\n");
283:     DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);
284:     DMTSCopy(oldsdm,sdm);
285:     DMTSDestroy((DMTS*)&dm->dmts);
286:     dm->dmts = (PetscObject) sdm;
287:   }
288:   *tsdm = sdm;
289:   return(0);
290: }

294: /*@C
295:    DMCopyDMTS - copies a DM context to a new DM

297:    Logically Collective

299:    Input Arguments:
300: +  dmsrc - DM to obtain context from
301: -  dmdest - DM to add context to

303:    Level: developer

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

308: .seealso: DMGetDMTS(), TSSetDM()
309: @*/
310: PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
311: {

317:   DMTSDestroy((DMTS*)&dmdest->dmts);
318:   dmdest->dmts = dmsrc->dmts;
319:   PetscObjectReference(dmdest->dmts);
320:   DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);
321:   DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);
322:   return(0);
323: }

327: /*@C
328:    DMTSSetIFunction - set TS implicit function evaluation function

330:    Not Collective

332:    Input Arguments:
333: +  dm - DM to be used with TS
334: .  func - function evaluation function, see TSSetIFunction() for calling sequence
335: -  ctx - context for residual evaluation

337:    Level: advanced

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

344: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
345: @*/
346: PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
347: {
349:   DMTS           tsdm;

353:   DMGetDMTSWrite(dm,&tsdm);
354:   if (func) tsdm->ops->ifunction = func;
355:   if (ctx)  tsdm->ifunctionctx = ctx;
356:   return(0);
357: }

361: /*@C
362:    DMTSGetIFunction - get TS implicit residual evaluation function

364:    Not Collective

366:    Input Argument:
367: .  dm - DM to be used with TS

369:    Output Arguments:
370: +  func - function evaluation function, see TSSetIFunction() for calling sequence
371: -  ctx - context for residual evaluation

373:    Level: advanced

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

379: .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
380: @*/
381: PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
382: {
384:   DMTS           tsdm;

388:   DMGetDMTS(dm,&tsdm);
389:   if (func) *func = tsdm->ops->ifunction;
390:   if (ctx)  *ctx = tsdm->ifunctionctx;
391:   return(0);
392: }

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

399:    Not Collective

401:    Input Arguments:
402: +  dm - DM to be used with TS
403: .  fun - function evaluation function, see TSSetI2Function() for calling sequence
404: -  ctx - context for residual evaluation

406:    Level: advanced

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

412: .seealso: TSSetI2Function()
413: @*/
414: PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx)
415: {
416:   DMTS           tsdm;

421:   DMGetDMTSWrite(dm,&tsdm);
422:   if (fun) tsdm->ops->i2function = fun;
423:   if (ctx) tsdm->i2functionctx   = ctx;
424:   return(0);
425: }

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

432:    Not Collective

434:    Input Argument:
435: .  dm - DM to be used with TS

437:    Output Arguments:
438: +  fun - function evaluation function, see TSSetI2Function() for calling sequence
439: -  ctx - context for residual evaluation

441:    Level: advanced

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

447: .seealso: DMTSSetI2Function(),TSGetI2Function()
448: @*/
449: PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx)
450: {
451:   DMTS           tsdm;

456:   DMGetDMTS(dm,&tsdm);
457:   if (fun) *fun = tsdm->ops->i2function;
458:   if (ctx) *ctx = tsdm->i2functionctx;
459:   return(0);
460: }

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

467:    Not Collective

469:    Input Arguments:
470: +  dm - DM to be used with TS
471: .  fun - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
472: -  ctx - context for Jacobian evaluation

474:    Level: advanced

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

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

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

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

500:    Not Collective

502:    Input Argument:
503: .  dm - DM to be used with TS

505:    Output Arguments:
506: +  jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
507: -  ctx - context for Jacobian evaluation

509:    Level: advanced

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

515: .seealso: DMTSSetI2Jacobian(),TSGetI2Jacobian()
516: @*/
517: PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
518: {
519:   DMTS           tsdm;

524:   DMGetDMTS(dm,&tsdm);
525:   if (jac) *jac = tsdm->ops->i2jacobian;
526:   if (ctx) *ctx = tsdm->i2jacobianctx;
527:   return(0);
528: }

532: /*@C
533:    DMTSSetRHSFunction - set TS explicit residual evaluation function

535:    Not Collective

537:    Input Arguments:
538: +  dm - DM to be used with TS
539: .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
540: -  ctx - context for residual evaluation

542:    Level: advanced

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

549: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
550: @*/
551: PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
552: {
554:   DMTS           tsdm;

558:   DMGetDMTSWrite(dm,&tsdm);
559:   if (func) tsdm->ops->rhsfunction = func;
560:   if (ctx)  tsdm->rhsfunctionctx = ctx;
561:   return(0);
562: }

566: /*@C
567:    DMTSGetSolutionFunction - gets the TS solution evaluation function

569:    Not Collective

571:    Input Arguments:
572: .  dm - DM to be used with TS

574:    Output Parameters:
575: +  func - solution function evaluation function, see TSSetSolution() for calling sequence
576: -  ctx - context for solution evaluation

578:    Level: advanced

580: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
581: @*/
582: PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
583: {
585:   DMTS           tsdm;

589:   DMGetDMTS(dm,&tsdm);
590:   if (func) *func = tsdm->ops->solution;
591:   if (ctx)  *ctx  = tsdm->solutionctx;
592:   return(0);
593: }

597: /*@C
598:    DMTSSetSolutionFunction - set TS solution evaluation function

600:    Not Collective

602:    Input Arguments:
603: +  dm - DM to be used with TS
604: .  func - solution function evaluation function, see TSSetSolution() for calling sequence
605: -  ctx - context for solution evaluation

607:    Level: advanced

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

614: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
615: @*/
616: PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
617: {
619:   DMTS           tsdm;

623:   DMGetDMTSWrite(dm,&tsdm);
624:   if (func) tsdm->ops->solution = func;
625:   if (ctx)  tsdm->solutionctx   = ctx;
626:   return(0);
627: }

631: /*@C
632:    DMTSSetForcingFunction - set TS forcing function evaluation function

634:    Not Collective

636:    Input Arguments:
637: +  dm - DM to be used with TS
638: .  f - forcing function evaluation function; see TSForcingFunction
639: -  ctx - context for solution evaluation

641:    Level: advanced

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

648: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
649: @*/
650: PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
651: {
653:   DMTS           tsdm;

657:   DMGetDMTSWrite(dm,&tsdm);
658:   if (f)    tsdm->ops->forcing = f;
659:   if (ctx)  tsdm->forcingctx   = ctx;
660:   return(0);
661: }


666: /*@C
667:    DMTSGetForcingFunction - get TS forcing function evaluation function

669:    Not Collective

671:    Input Argument:
672: .   dm - DM to be used with TS

674:    Output Arguments:
675: +  f - forcing function evaluation function; see TSForcingFunction for details
676: -  ctx - context for solution evaluation

678:    Level: advanced

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

685: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
686: @*/
687: PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
688: {
690:   DMTS           tsdm;

694:   DMGetDMTSWrite(dm,&tsdm);
695:   if (f)   *f   = tsdm->ops->forcing;
696:   if (ctx) *ctx = tsdm->forcingctx;
697:   return(0);
698: }

702: /*@C
703:    DMTSGetRHSFunction - get TS explicit residual evaluation function

705:    Not Collective

707:    Input Argument:
708: .  dm - DM to be used with TS

710:    Output Arguments:
711: +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
712: -  ctx - context for residual evaluation

714:    Level: advanced

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

720: .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
721: @*/
722: PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
723: {
725:   DMTS           tsdm;

729:   DMGetDMTS(dm,&tsdm);
730:   if (func) *func = tsdm->ops->rhsfunction;
731:   if (ctx)  *ctx = tsdm->rhsfunctionctx;
732:   return(0);
733: }

737: /*@C
738:    DMTSSetIJacobian - set TS Jacobian evaluation function

740:    Not Collective

742:    Input Argument:
743: +  dm - DM to be used with TS
744: .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
745: -  ctx - context for residual evaluation

747:    Level: advanced

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

754: .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
755: @*/
756: PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
757: {
759:   DMTS           sdm;

763:   DMGetDMTSWrite(dm,&sdm);
764:   if (func) sdm->ops->ijacobian = func;
765:   if (ctx)  sdm->ijacobianctx   = ctx;
766:   return(0);
767: }

771: /*@C
772:    DMTSGetIJacobian - get TS Jacobian evaluation function

774:    Not Collective

776:    Input Argument:
777: .  dm - DM to be used with TS

779:    Output Arguments:
780: +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
781: -  ctx - context for residual evaluation

783:    Level: advanced

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

790: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
791: @*/
792: PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
793: {
795:   DMTS           tsdm;

799:   DMGetDMTS(dm,&tsdm);
800:   if (func) *func = tsdm->ops->ijacobian;
801:   if (ctx)  *ctx = tsdm->ijacobianctx;
802:   return(0);
803: }


808: /*@C
809:    DMTSSetRHSJacobian - set TS Jacobian evaluation function

811:    Not Collective

813:    Input Argument:
814: +  dm - DM to be used with TS
815: .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
816: -  ctx - context for residual evaluation

818:    Level: advanced

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

825: .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
826: @*/
827: PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
828: {
830:   DMTS           tsdm;

834:   DMGetDMTSWrite(dm,&tsdm);
835:   if (func) tsdm->ops->rhsjacobian = func;
836:   if (ctx)  tsdm->rhsjacobianctx = ctx;
837:   return(0);
838: }

842: /*@C
843:    DMTSGetRHSJacobian - get TS Jacobian evaluation function

845:    Not Collective

847:    Input Argument:
848: .  dm - DM to be used with TS

850:    Output Arguments:
851: +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
852: -  ctx - context for residual evaluation

854:    Level: advanced

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

861: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
862: @*/
863: PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
864: {
866:   DMTS           tsdm;

870:   DMGetDMTS(dm,&tsdm);
871:   if (func) *func = tsdm->ops->rhsjacobian;
872:   if (ctx)  *ctx = tsdm->rhsjacobianctx;
873:   return(0);
874: }

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

881:    Not Collective

883:    Input Arguments:
884: +  dm - DM to be used with TS
885: .  view - viewer function
886: -  load - loading function

888:    Level: advanced

890: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
891: @*/
892: PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
893: {
895:   DMTS           tsdm;

899:   DMGetDMTSWrite(dm,&tsdm);
900:   tsdm->ops->ifunctionview = view;
901:   tsdm->ops->ifunctionload = load;
902:   return(0);
903: }

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

910:    Not Collective

912:    Input Arguments:
913: +  dm - DM to be used with TS
914: .  view - viewer function
915: -  load - loading function

917:    Level: advanced

919: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
920: @*/
921: PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
922: {
924:   DMTS           tsdm;

928:   DMGetDMTSWrite(dm,&tsdm);
929:   tsdm->ops->ijacobianview = view;
930:   tsdm->ops->ijacobianload = load;
931:   return(0);
932: }