Actual source code: dmts.c

petsc-3.4.5 2014-06-29
  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,PETSC_FUNCTION);
 27:   PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,PETSC_FUNCTION);
 28:   PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,PETSC_FUNCTION);
 29:   if (kdm->ops->ifunctionload) {
 30:     (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);
 31:   }
 32:   PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,PETSC_FUNCTION);
 33:   PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,PETSC_FUNCTION);
 34:   PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,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:       PetscErrorCode (*ifunctionview)(void*,PetscViewer);
 68:       PetscErrorCode (*ifunctionload)(void**,PetscViewer);
 69:     } funcstruct = {kdm->ops->ifunction,
 70:                     kdm->ops->ifunctionview,
 71:                     kdm->ops->ifunctionload};
 72:     struct {
 73:       TSIJacobian ijacobian;
 74:       PetscErrorCode (*ijacobianview)(void*,PetscViewer);
 75:       PetscErrorCode (*ijacobianload)(void**,PetscViewer);
 76:     } jacstruct = {kdm->ops->ijacobian,
 77:                    kdm->ops->ijacobianview,
 78:                    kdm->ops->ijacobianload};

 80:     PetscViewerBinaryWrite(viewer,&funcstruct,3,PETSC_FUNCTION,PETSC_FALSE);
 81:     if (kdm->ops->ifunctionview) {
 82:       (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);
 83:     }
 84:     PetscViewerBinaryWrite(viewer,&jacstruct,3,PETSC_FUNCTION,PETSC_FALSE);
 85:     if (kdm->ops->ijacobianview) {
 86:       (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);
 87:     }
 88:   }
 89:   return(0);
 90: }

 94: static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
 95: {

 99: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
100:   TSInitializePackage();
101: #endif
102:   PetscHeaderCreate(*kdm, _p_DMTS, struct _DMTSOps, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);
103:   PetscMemzero((*kdm)->ops, sizeof(struct _DMTSOps));
104:   return(0);
105: }

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

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

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

129:   return(0);
130: }

134: static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
135: {

139:   DMCopyDMTS(dm,subdm);
140:   return(0);
141: }

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

155: /*@C
156:    DMTSCopy - copies the information in a DMTS to another DMTS

158:    Not Collective

160:    Input Argument:
161: +  kdm - Original DMTS
162: -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()

164:    Level: developer

166: .seealso: DMTSCreate(), DMTSDestroy()
167: @*/
168: PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
169: {

175:   nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
176:   nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
177:   nkdm->ops->ifunction   = kdm->ops->ifunction;
178:   nkdm->ops->ijacobian   = kdm->ops->ijacobian;
179:   nkdm->ops->solution    = kdm->ops->solution;
180:   nkdm->ops->destroy     = kdm->ops->destroy;
181:   nkdm->ops->duplicate   = kdm->ops->duplicate;

183:   nkdm->rhsfunctionctx = kdm->rhsfunctionctx;
184:   nkdm->rhsjacobianctx = kdm->rhsjacobianctx;
185:   nkdm->ifunctionctx   = kdm->ifunctionctx;
186:   nkdm->ijacobianctx   = kdm->ijacobianctx;
187:   nkdm->solutionctx    = kdm->solutionctx;

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

191:   /*
192:   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
193:   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
194:   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
195:   */

197:   /* implementation specific copy hooks */
198:   if (kdm->ops->duplicate) {(*kdm->ops->duplicate)(kdm,nkdm);}
199:   return(0);
200: }

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

207:    Not Collective

209:    Input Argument:
210: .  dm - DM to be used with TS

212:    Output Argument:
213: .  tsdm - private DMTS context

215:    Level: developer

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

220: .seealso: DMGetDMTSWrite()
221: @*/
222: PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
223: {

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

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

244:    Not Collective

246:    Input Argument:
247: .  dm - DM to be used with TS

249:    Output Argument:
250: .  tsdm - private DMTS context

252:    Level: developer

254: .seealso: DMGetDMTS()
255: @*/
256: PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
257: {
259:   DMTS           sdm;

263:   DMGetDMTS(dm,&sdm);
264:   if (!sdm->originaldm) sdm->originaldm = dm;
265:   if (sdm->originaldm != dm) {  /* Copy on write */
266:     DMTS oldsdm = sdm;
267:     PetscInfo(dm,"Copying DMTS due to write\n");
268:     DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);
269:     DMTSCopy(oldsdm,sdm);
270:     DMTSDestroy((DMTS*)&dm->dmts);
271:     dm->dmts = (PetscObject) sdm;
272:   }
273:   *tsdm = sdm;
274:   return(0);
275: }

279: /*@C
280:    DMCopyDMTS - copies a DM context to a new DM

282:    Logically Collective

284:    Input Arguments:
285: +  dmsrc - DM to obtain context from
286: -  dmdest - DM to add context to

288:    Level: developer

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

293: .seealso: DMGetDMTS(), TSSetDM()
294: @*/
295: PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
296: {

302:   DMTSDestroy((DMTS*)&dmdest->dmts);
303:   dmdest->dmts = dmsrc->dmts;
304:   PetscObjectReference(dmdest->dmts);
305:   DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);
306:   DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);
307:   return(0);
308: }

312: /*@C
313:    DMTSSetIFunction - set TS implicit function evaluation function

315:    Not Collective

317:    Input Arguments:
318: +  dm - DM to be used with TS
319: .  func - function evaluation function, see TSSetIFunction() for calling sequence
320: -  ctx - context for residual evaluation

322:    Level: advanced

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.

329: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
330: @*/
331: PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
332: {
334:   DMTS           tsdm;

338:   DMGetDMTSWrite(dm,&tsdm);
339:   if (func) tsdm->ops->ifunction = func;
340:   if (ctx)  tsdm->ifunctionctx = ctx;
341:   return(0);
342: }

346: /*@C
347:    DMTSGetIFunction - get TS implicit residual evaluation function

349:    Not Collective

351:    Input Argument:
352: .  dm - DM to be used with TS

354:    Output Arguments:
355: +  func - function evaluation function, see TSSetIFunction() for calling sequence
356: -  ctx - context for residual evaluation

358:    Level: advanced

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

364: .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
365: @*/
366: PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
367: {
369:   DMTS           tsdm;

373:   DMGetDMTS(dm,&tsdm);
374:   if (func) *func = tsdm->ops->ifunction;
375:   if (ctx)  *ctx = tsdm->ifunctionctx;
376:   return(0);
377: }


382: /*@C
383:    DMTSSetRHSFunction - set TS explicit residual evaluation function

385:    Not Collective

387:    Input Arguments:
388: +  dm - DM to be used with TS
389: .  func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence
390: -  ctx - context for residual evaluation

392:    Level: advanced

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

399: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
400: @*/
401: PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
402: {
404:   DMTS           tsdm;

408:   DMGetDMTSWrite(dm,&tsdm);
409:   if (func) tsdm->ops->rhsfunction = func;
410:   if (ctx)  tsdm->rhsfunctionctx = ctx;
411:   return(0);
412: }

416: /*@C
417:    DMTSGetSolutionFunction - gets the TS solution evaluation function

419:    Not Collective

421:    Input Arguments:
422: .  dm - DM to be used with TS

424:    Output Parameters:
425: +  func - solution function evaluation function, see TSSetSolution() for calling sequence
426: -  ctx - context for solution evaluation

428:    Level: advanced

430: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
431: @*/
432: PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
433: {
435:   DMTS           tsdm;

439:   DMGetDMTS(dm,&tsdm);
440:   if (func) *func = tsdm->ops->solution;
441:   if (ctx)  *ctx  = tsdm->solutionctx;
442:   return(0);
443: }

447: /*@C
448:    DMTSSetSolutionFunction - set TS solution evaluation function

450:    Not Collective

452:    Input Arguments:
453: +  dm - DM to be used with TS
454: .  func - solution function evaluation function, see TSSetSolution() for calling sequence
455: -  ctx - context for solution evaluation

457:    Level: advanced

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

464: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
465: @*/
466: PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
467: {
469:   DMTS           tsdm;

473:   DMGetDMTSWrite(dm,&tsdm);
474:   if (func) tsdm->ops->solution = func;
475:   if (ctx)  tsdm->solutionctx   = ctx;
476:   return(0);
477: }

481: /*@C
482:    DMTSSetForcingFunction - set TS forcing function evaluation function

484:    Not Collective

486:    Input Arguments:
487: +  dm - DM to be used with TS
488: .  TSForcingFunction - forcing function evaluation function
489: -  ctx - context for solution evaluation

491:    Level: advanced

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

498: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
499: @*/
500: PetscErrorCode DMTSSetForcingFunction(DM dm,PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*),void *ctx)
501: {
503:   DMTS           tsdm;

507:   DMGetDMTSWrite(dm,&tsdm);
508:   if (TSForcingFunction) tsdm->ops->forcing = TSForcingFunction;
509:   if (ctx)  tsdm->forcingctx   = ctx;
510:   return(0);
511: }


516: /*@C
517:    DMTSGetForcingFunction - get TS forcing function evaluation function

519:    Not Collective

521:    Input Argument:
522: .   dm - DM to be used with TS

524:    Output Arguments:
525: +  TSForcingFunction - forcing function evaluation function
526: -  ctx - context for solution evaluation

528:    Level: advanced

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

535: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
536: @*/
537: PetscErrorCode DMTSGetForcingFunction(DM dm,PetscErrorCode (**TSForcingFunction)(TS,PetscReal,Vec,void*),void **ctx)
538: {
540:   DMTS           tsdm;

544:   DMGetDMTSWrite(dm,&tsdm);
545:   if (TSForcingFunction) *TSForcingFunction = tsdm->ops->forcing;
546:   if (ctx) *ctx = tsdm->forcingctx;
547:   return(0);
548: }

552: /*@C
553:    DMTSGetRHSFunction - get TS explicit residual evaluation function

555:    Not Collective

557:    Input Argument:
558: .  dm - DM to be used with TS

560:    Output Arguments:
561: +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
562: -  ctx - context for residual evaluation

564:    Level: advanced

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

570: .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
571: @*/
572: PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
573: {
575:   DMTS           tsdm;

579:   DMGetDMTS(dm,&tsdm);
580:   if (func) *func = tsdm->ops->rhsfunction;
581:   if (ctx)  *ctx = tsdm->rhsfunctionctx;
582:   return(0);
583: }

587: /*@C
588:    DMTSSetIJacobian - set TS Jacobian evaluation function

590:    Not Collective

592:    Input Argument:
593: +  dm - DM to be used with TS
594: .  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
595: -  ctx - context for residual evaluation

597:    Level: advanced

599:    Note:
600:    TSSetJacobian() 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 Jacobian.

604: .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
605: @*/
606: PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
607: {
609:   DMTS           sdm;

613:   DMGetDMTSWrite(dm,&sdm);
614:   if (func) sdm->ops->ijacobian = func;
615:   if (ctx)  sdm->ijacobianctx   = ctx;
616:   return(0);
617: }

621: /*@C
622:    DMTSGetIJacobian - get TS Jacobian evaluation function

624:    Not Collective

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

629:    Output Arguments:
630: +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
631: -  ctx - context for residual evaluation

633:    Level: advanced

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

640: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
641: @*/
642: PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
643: {
645:   DMTS           tsdm;

649:   DMGetDMTS(dm,&tsdm);
650:   if (func) *func = tsdm->ops->ijacobian;
651:   if (ctx)  *ctx = tsdm->ijacobianctx;
652:   return(0);
653: }


658: /*@C
659:    DMTSSetRHSJacobian - set TS Jacobian evaluation function

661:    Not Collective

663:    Input Argument:
664: +  dm - DM to be used with TS
665: .  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
666: -  ctx - context for residual evaluation

668:    Level: advanced

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

675: .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian()
676: @*/
677: PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
678: {
680:   DMTS           tsdm;

684:   DMGetDMTSWrite(dm,&tsdm);
685:   if (func) tsdm->ops->rhsjacobian = func;
686:   if (ctx)  tsdm->rhsjacobianctx = ctx;
687:   return(0);
688: }

692: /*@C
693:    DMTSGetRHSJacobian - get TS Jacobian evaluation function

695:    Not Collective

697:    Input Argument:
698: .  dm - DM to be used with TS

700:    Output Arguments:
701: +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
702: -  ctx - context for residual evaluation

704:    Level: advanced

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

711: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
712: @*/
713: PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
714: {
716:   DMTS           tsdm;

720:   DMGetDMTS(dm,&tsdm);
721:   if (func) *func = tsdm->ops->rhsjacobian;
722:   if (ctx)  *ctx = tsdm->rhsjacobianctx;
723:   return(0);
724: }

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

731:    Not Collective

733:    Input Arguments:
734: +  dm - DM to be used with TS
735: .  view - viewer function
736: -  load - loading function

738:    Level: advanced

740: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
741: @*/
742: PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
743: {
745:   DMTS           tsdm;

749:   DMGetDMTSWrite(dm,&tsdm);
750:   tsdm->ops->ifunctionview = view;
751:   tsdm->ops->ifunctionload = load;
752:   return(0);
753: }

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

760:    Not Collective

762:    Input Arguments:
763: +  dm - DM to be used with TS
764: .  view - viewer function
765: -  load - loading function

767:    Level: advanced

769: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
770: @*/
771: PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
772: {
774:   DMTS           tsdm;

778:   DMGetDMTSWrite(dm,&tsdm);
779:   tsdm->ops->ijacobianview = view;
780:   tsdm->ops->ijacobianload = load;
781:   return(0);
782: }