Actual source code: dmts.c

petsc-3.9.4 2018-09-11
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,PETSC_FALSE);
 82:     PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION,PETSC_FALSE);
 83:     PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION,PETSC_FALSE);
 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,PETSC_FALSE);
 91:     PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION,PETSC_FALSE);
 92:     PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION,PETSC_FALSE);
 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:     DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);
229:     DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);
230:   }
231:   return(0);
232: }

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

237:    Not Collective

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

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

245:    Level: developer

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

256:   DMGetDMTS(dm,&sdm);
257:   if (!sdm->originaldm) sdm->originaldm = dm;
258:   if (sdm->originaldm != dm) {  /* Copy on write */
259:     DMTS oldsdm = sdm;
260:     PetscInfo(dm,"Copying DMTS due to write\n");
261:     DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);
262:     DMTSCopy(oldsdm,sdm);
263:     DMTSDestroy((DMTS*)&dm->dmts);
264:     dm->dmts = (PetscObject) sdm;
265:   }
266:   *tsdm = sdm;
267:   return(0);
268: }

270: /*@C
271:    DMCopyDMTS - copies a DM context to a new DM

273:    Logically Collective

275:    Input Arguments:
276: +  dmsrc - DM to obtain context from
277: -  dmdest - DM to add context to

279:    Level: developer

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

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

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

301: /*@C
302:    DMTSSetIFunction - set TS implicit function evaluation function

304:    Not Collective

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

311:    Level: advanced

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

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

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

333: /*@C
334:    DMTSGetIFunction - get TS implicit residual evaluation function

336:    Not Collective

338:    Input Argument:
339: .  dm - DM to be used with TS

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

345:    Level: advanced

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

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

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

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

369:    Not Collective

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

376:    Level: advanced

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

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

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

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

400:    Not Collective

402:    Input Argument:
403: .  dm - DM to be used with TS

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

409:    Level: advanced

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

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

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

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

433:    Not Collective

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

440:    Level: advanced

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

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

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

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

464:    Not Collective

466:    Input Argument:
467: .  dm - DM to be used with TS

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

473:    Level: advanced

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

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

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

494: /*@C
495:    DMTSSetRHSFunction - set TS explicit residual evaluation function

497:    Not Collective

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

504:    Level: advanced

506:    Note:
507:    TSSetRSHFunction() 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.

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

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

526: /*@C
527:    DMTSGetSolutionFunction - gets the TS solution evaluation function

529:    Not Collective

531:    Input Arguments:
532: .  dm - DM to be used with TS

534:    Output Parameters:
535: +  func - solution function evaluation function, see TSSetSolution() for calling sequence
536: -  ctx - context for solution evaluation

538:    Level: advanced

540: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSSetSolutionFunction()
541: @*/
542: PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
543: {
545:   DMTS           tsdm;

549:   DMGetDMTS(dm,&tsdm);
550:   if (func) *func = tsdm->ops->solution;
551:   if (ctx)  *ctx  = tsdm->solutionctx;
552:   return(0);
553: }

555: /*@C
556:    DMTSSetSolutionFunction - set TS solution evaluation function

558:    Not Collective

560:    Input Arguments:
561: +  dm - DM to be used with TS
562: .  func - solution function evaluation function, see TSSetSolution() for calling sequence
563: -  ctx - context for solution evaluation

565:    Level: advanced

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

572: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSGetSolutionFunction()
573: @*/
574: PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
575: {
577:   DMTS           tsdm;

581:   DMGetDMTSWrite(dm,&tsdm);
582:   if (func) tsdm->ops->solution = func;
583:   if (ctx)  tsdm->solutionctx   = ctx;
584:   return(0);
585: }

587: /*@C
588:    DMTSSetForcingFunction - set TS forcing function evaluation function

590:    Not Collective

592:    Input Arguments:
593: +  dm - DM to be used with TS
594: .  f - forcing function evaluation function; see TSForcingFunction
595: -  ctx - context for solution evaluation

597:    Level: advanced

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

604: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
605: @*/
606: PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
607: {
609:   DMTS           tsdm;

613:   DMGetDMTSWrite(dm,&tsdm);
614:   if (f)    tsdm->ops->forcing = f;
615:   if (ctx)  tsdm->forcingctx   = ctx;
616:   return(0);
617: }


620: /*@C
621:    DMTSGetForcingFunction - get TS forcing function evaluation function

623:    Not Collective

625:    Input Argument:
626: .   dm - DM to be used with TS

628:    Output Arguments:
629: +  f - forcing function evaluation function; see TSForcingFunction for details
630: -  ctx - context for solution evaluation

632:    Level: advanced

634:    Note:
635:    TSSetForcingFunction() 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(), TSSetForcingFunction(), DMTSGetForcingFunction()
640: @*/
641: PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
642: {
644:   DMTS           tsdm;

648:   DMGetDMTSWrite(dm,&tsdm);
649:   if (f)   *f   = tsdm->ops->forcing;
650:   if (ctx) *ctx = tsdm->forcingctx;
651:   return(0);
652: }

654: /*@C
655:    DMTSGetRHSFunction - get TS explicit residual evaluation function

657:    Not Collective

659:    Input Argument:
660: .  dm - DM to be used with TS

662:    Output Arguments:
663: +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
664: -  ctx - context for residual evaluation

666:    Level: advanced

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

672: .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
673: @*/
674: PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
675: {
677:   DMTS           tsdm;

681:   DMGetDMTS(dm,&tsdm);
682:   if (func) *func = tsdm->ops->rhsfunction;
683:   if (ctx)  *ctx = tsdm->rhsfunctionctx;
684:   return(0);
685: }

687: /*@C
688:    DMTSSetIJacobian - set TS Jacobian evaluation function

690:    Not Collective

692:    Input Argument:
693: +  dm - DM to be used with TS
694: .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
695: -  ctx - context for residual evaluation

697:    Level: advanced

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

704: .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
705: @*/
706: PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
707: {
709:   DMTS           sdm;

713:   DMGetDMTSWrite(dm,&sdm);
714:   if (func) sdm->ops->ijacobian = func;
715:   if (ctx)  sdm->ijacobianctx   = ctx;
716:   return(0);
717: }

719: /*@C
720:    DMTSGetIJacobian - get TS Jacobian evaluation function

722:    Not Collective

724:    Input Argument:
725: .  dm - DM to be used with TS

727:    Output Arguments:
728: +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
729: -  ctx - context for residual evaluation

731:    Level: advanced

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

738: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
739: @*/
740: PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
741: {
743:   DMTS           tsdm;

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


754: /*@C
755:    DMTSSetRHSJacobian - 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 TSSetRHSJacobian() 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 DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
774: {
776:   DMTS           tsdm;

780:   DMGetDMTSWrite(dm,&tsdm);
781:   if (func) tsdm->ops->rhsjacobian = func;
782:   if (ctx)  tsdm->rhsjacobianctx = ctx;
783:   return(0);
784: }

786: /*@C
787:    DMTSGetRHSJacobian - 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 TSSetRHSJacobian() 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 DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
808: {
810:   DMTS           tsdm;

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

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

823:    Not Collective

825:    Input Arguments:
826: +  dm - DM to be used with TS
827: .  view - viewer function
828: -  load - loading function

830:    Level: advanced

832: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
833: @*/
834: PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
835: {
837:   DMTS           tsdm;

841:   DMGetDMTSWrite(dm,&tsdm);
842:   tsdm->ops->ifunctionview = view;
843:   tsdm->ops->ifunctionload = load;
844:   return(0);
845: }

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

850:    Not Collective

852:    Input Arguments:
853: +  dm - DM to be used with TS
854: .  view - viewer function
855: -  load - loading function

857:    Level: advanced

859: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
860: @*/
861: PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
862: {
864:   DMTS           tsdm;

868:   DMGetDMTSWrite(dm,&tsdm);
869:   tsdm->ops->ijacobianview = view;
870:   tsdm->ops->ijacobianload = load;
871:   return(0);
872: }