Actual source code: dmts.c

  1: #include <petsc/private/tsimpl.h>
  2: #include <petsc/private/dmimpl.h>

  4: static PetscErrorCode DMTSDestroy(DMTS *kdm)
  5: {
  6:   if (!*kdm) return 0;
  8:   if (--((PetscObject)(*kdm))->refct > 0) {*kdm = NULL; return 0;}
  9:   if ((*kdm)->ops->destroy) ((*kdm)->ops->destroy)(*kdm);
 10:   PetscHeaderDestroy(kdm);
 11:   return 0;
 12: }

 14: PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer)
 15: {
 16:   PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,NULL,PETSC_FUNCTION);
 17:   PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,NULL,PETSC_FUNCTION);
 18:   PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,NULL,PETSC_FUNCTION);
 19:   if (kdm->ops->ifunctionload) {
 20:     (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);
 21:   }
 22:   PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,NULL,PETSC_FUNCTION);
 23:   PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,NULL,PETSC_FUNCTION);
 24:   PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,NULL,PETSC_FUNCTION);
 25:   if (kdm->ops->ijacobianload) {
 26:     (*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer);
 27:   }
 28:   return 0;
 29: }

 31: PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
 32: {
 33:   PetscBool      isascii,isbinary;

 35:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 36:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
 37:   if (isascii) {
 38: #if defined(PETSC_SERIALIZE_FUNCTIONS)
 39:     const char *fname;

 41:     PetscFPTFind(kdm->ops->ifunction,&fname);
 42:     if (fname) {
 43:       PetscViewerASCIIPrintf(viewer,"  IFunction used by TS: %s\n",fname);
 44:     }
 45:     PetscFPTFind(kdm->ops->ijacobian,&fname);
 46:     if (fname) {
 47:       PetscViewerASCIIPrintf(viewer,"  IJacobian function used by TS: %s\n",fname);
 48:     }
 49: #endif
 50:   } else if (isbinary) {
 51:     struct {
 52:       TSIFunction ifunction;
 53:     } funcstruct;
 54:     struct {
 55:       PetscErrorCode (*ifunctionview)(void*,PetscViewer);
 56:     } funcviewstruct;
 57:     struct {
 58:       PetscErrorCode (*ifunctionload)(void**,PetscViewer);
 59:     } funcloadstruct;
 60:     struct {
 61:       TSIJacobian ijacobian;
 62:     } jacstruct;
 63:     struct {
 64:       PetscErrorCode (*ijacobianview)(void*,PetscViewer);
 65:     } jacviewstruct;
 66:     struct {
 67:       PetscErrorCode (*ijacobianload)(void**,PetscViewer);
 68:     } jacloadstruct;

 70:     funcstruct.ifunction         = kdm->ops->ifunction;
 71:     funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
 72:     funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
 73:     PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION);
 74:     PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION);
 75:     PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION);
 76:     if (kdm->ops->ifunctionview) {
 77:       (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);
 78:     }
 79:     jacstruct.ijacobian = kdm->ops->ijacobian;
 80:     jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
 81:     jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
 82:     PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION);
 83:     PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION);
 84:     PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION);
 85:     if (kdm->ops->ijacobianview) {
 86:       (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);
 87:     }
 88:   }
 89:   return 0;
 90: }

 92: static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
 93: {
 94:   TSInitializePackage();
 95:   PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);
 96:   return 0;
 97: }

 99: /* Attaches the DMTS to the coarse level.
100:  * Under what conditions should we copy versus duplicate?
101:  */
102: static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
103: {
104:   DMCopyDMTS(dm,dmc);
105:   return 0;
106: }

108: /* This could restrict auxiliary information to the coarse level.
109:  */
110: static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
111: {
112:   return 0;
113: }

115: static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
116: {
117:   DMCopyDMTS(dm,subdm);
118:   return 0;
119: }

121: /* This could restrict auxiliary information to the coarse level.
122:  */
123: static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
124: {
125:   return 0;
126: }

128: /*@C
129:    DMTSCopy - copies the information in a DMTS to another DMTS

131:    Not Collective

133:    Input Parameters:
134: +  kdm - Original DMTS
135: -  nkdm - DMTS to receive the data, should have been created with DMTSCreate()

137:    Level: developer

139: .seealso: DMTSCreate(), DMTSDestroy()
140: @*/
141: PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
142: {
145:   nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
146:   nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
147:   nkdm->ops->ifunction   = kdm->ops->ifunction;
148:   nkdm->ops->ijacobian   = kdm->ops->ijacobian;
149:   nkdm->ops->i2function  = kdm->ops->i2function;
150:   nkdm->ops->i2jacobian  = kdm->ops->i2jacobian;
151:   nkdm->ops->solution    = kdm->ops->solution;
152:   nkdm->ops->destroy     = kdm->ops->destroy;
153:   nkdm->ops->duplicate   = kdm->ops->duplicate;

155:   nkdm->rhsfunctionctx = kdm->rhsfunctionctx;
156:   nkdm->rhsjacobianctx = kdm->rhsjacobianctx;
157:   nkdm->ifunctionctx   = kdm->ifunctionctx;
158:   nkdm->ijacobianctx   = kdm->ijacobianctx;
159:   nkdm->i2functionctx  = kdm->i2functionctx;
160:   nkdm->i2jacobianctx  = kdm->i2jacobianctx;
161:   nkdm->solutionctx    = kdm->solutionctx;

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

165:   /*
166:   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
167:   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
168:   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
169:   */

171:   /* implementation specific copy hooks */
172:   if (kdm->ops->duplicate) (*kdm->ops->duplicate)(kdm,nkdm);
173:   return 0;
174: }

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

179:    Not Collective

181:    Input Parameter:
182: .  dm - DM to be used with TS

184:    Output Parameter:
185: .  tsdm - private DMTS context

187:    Level: developer

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

192: .seealso: DMGetDMTSWrite()
193: @*/
194: PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
195: {
197:   *tsdm = (DMTS) dm->dmts;
198:   if (!*tsdm) {
199:     PetscInfo(dm,"Creating new DMTS\n");
200:     DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm);
201:     dm->dmts = (PetscObject) *tsdm;
202:     (*tsdm)->originaldm = dm;
203:     DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);
204:     DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);
205:   }
206:   return 0;
207: }

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

212:    Not Collective

214:    Input Parameter:
215: .  dm - DM to be used with TS

217:    Output Parameter:
218: .  tsdm - private DMTS context

220:    Level: developer

222: .seealso: DMGetDMTS()
223: @*/
224: PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
225: {
226:   DMTS           sdm;

229:   DMGetDMTS(dm,&sdm);
231:   if (sdm->originaldm != dm) {  /* Copy on write */
232:     DMTS oldsdm = sdm;
233:     PetscInfo(dm,"Copying DMTS due to write\n");
234:     DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);
235:     DMTSCopy(oldsdm,sdm);
236:     DMTSDestroy((DMTS*)&dm->dmts);
237:     dm->dmts = (PetscObject) sdm;
238:     sdm->originaldm = dm;
239:   }
240:   *tsdm = sdm;
241:   return 0;
242: }

244: /*@C
245:    DMCopyDMTS - copies a DM context to a new DM

247:    Logically Collective

249:    Input Parameters:
250: +  dmsrc - DM to obtain context from
251: -  dmdest - DM to add context to

253:    Level: developer

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

258: .seealso: DMGetDMTS(), TSSetDM()
259: @*/
260: PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
261: {
264:   DMTSDestroy((DMTS*)&dmdest->dmts);
265:   dmdest->dmts = dmsrc->dmts;
266:   PetscObjectReference(dmdest->dmts);
267:   DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);
268:   DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);
269:   return 0;
270: }

272: /*@C
273:    DMTSSetIFunction - set TS implicit function evaluation function

275:    Not Collective

277:    Input Parameters:
278: +  dm - DM to be used with TS
279: .  func - function evaluating f(t,u,u_t)
280: -  ctx - context for residual evaluation

282:    Calling sequence of func:
283: $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);

285: +  t   - time at step/stage being solved
286: .  u   - state vector
287: .  u_t - time derivative of state vector
288: .  F   - function vector
289: -  ctx - [optional] user-defined context for matrix evaluation routine

291:    Level: advanced

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

298: .seealso: DMTSSetContext(), TSSetIFunction(), DMTSSetJacobian()
299: @*/
300: PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx)
301: {
302:   DMTS           tsdm;

305:   DMGetDMTSWrite(dm,&tsdm);
306:   if (func) tsdm->ops->ifunction = func;
307:   if (ctx)  tsdm->ifunctionctx = ctx;
308:   return 0;
309: }

311: /*@C
312:    DMTSGetIFunction - get TS implicit residual evaluation function

314:    Not Collective

316:    Input Parameter:
317: .  dm - DM to be used with TS

319:    Output Parameters:
320: +  func - function evaluation function, see TSSetIFunction() for calling sequence
321: -  ctx - context for residual evaluation

323:    Level: advanced

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

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

336:   DMGetDMTS(dm,&tsdm);
337:   if (func) *func = tsdm->ops->ifunction;
338:   if (ctx)  *ctx = tsdm->ifunctionctx;
339:   return 0;
340: }

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

345:    Not Collective

347:    Input Parameters:
348: +  dm - DM to be used with TS
349: .  fun - function evaluation routine
350: -  ctx - context for residual evaluation

352:    Calling sequence of fun:
353: $     PetscErrorCode fun(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,Vec F,ctx);

355: +  t    - time at step/stage being solved
356: .  U    - state vector
357: .  U_t  - time derivative of state vector
358: .  U_tt - second time derivative of state vector
359: .  F    - function vector
360: -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)

362:    Level: advanced

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

368: .seealso: TSSetI2Function()
369: @*/
370: PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx)
371: {
372:   DMTS           tsdm;

375:   DMGetDMTSWrite(dm,&tsdm);
376:   if (fun) tsdm->ops->i2function = fun;
377:   if (ctx) tsdm->i2functionctx   = ctx;
378:   return 0;
379: }

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

384:    Not Collective

386:    Input Parameter:
387: .  dm - DM to be used with TS

389:    Output Parameters:
390: +  fun - function evaluation function, see TSSetI2Function() for calling sequence
391: -  ctx - context for residual evaluation

393:    Level: advanced

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

399: .seealso: DMTSSetI2Function(),TSGetI2Function()
400: @*/
401: PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx)
402: {
403:   DMTS           tsdm;

406:   DMGetDMTS(dm,&tsdm);
407:   if (fun) *fun = tsdm->ops->i2function;
408:   if (ctx) *ctx = tsdm->i2functionctx;
409:   return 0;
410: }

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

415:    Not Collective

417:    Input Parameters:
418: +  dm - DM to be used with TS
419: .  fun - Jacobian evaluation routine
420: -  ctx - context for Jacobian evaluation

422:    Calling sequence of jac:
423: $    PetscErrorCode jac(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,PetscReal v,PetscReal a,Mat J,Mat P,void *ctx);

425: +  t    - time at step/stage being solved
426: .  U    - state vector
427: .  U_t  - time derivative of state vector
428: .  U_tt - second time derivative of state vector
429: .  v    - shift for U_t
430: .  a    - shift for U_tt
431: .  J    - Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t  + a*dF/dU_tt
432: .  P    - preconditioning matrix for J, may be same as J
433: -  ctx  - [optional] user-defined context for matrix evaluation routine

435:    Level: advanced

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

441: .seealso: TSSetI2Jacobian()
442: @*/
443: PetscErrorCode DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void *ctx)
444: {
445:   DMTS           tsdm;

448:   DMGetDMTSWrite(dm,&tsdm);
449:   if (jac) tsdm->ops->i2jacobian = jac;
450:   if (ctx) tsdm->i2jacobianctx   = ctx;
451:   return 0;
452: }

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

457:    Not Collective

459:    Input Parameter:
460: .  dm - DM to be used with TS

462:    Output Parameters:
463: +  jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
464: -  ctx - context for Jacobian evaluation

466:    Level: advanced

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

472: .seealso: DMTSSetI2Jacobian(),TSGetI2Jacobian()
473: @*/
474: PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
475: {
476:   DMTS           tsdm;

479:   DMGetDMTS(dm,&tsdm);
480:   if (jac) *jac = tsdm->ops->i2jacobian;
481:   if (ctx) *ctx = tsdm->i2jacobianctx;
482:   return 0;
483: }

485: /*@C
486:    DMTSSetRHSFunction - set TS explicit residual evaluation function

488:    Not Collective

490:    Input Parameters:
491: +  dm - DM to be used with TS
492: .  func - RHS function evaluation routine
493: -  ctx - context for residual evaluation

495:     Calling sequence of func:
496: $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec F,void *ctx);

498: +   ts - timestep context
499: .   t - current timestep
500: .   u - input vector
501: .   F - function vector
502: -   ctx - [optional] user-defined function context

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(), TSSetRHSFunction(), DMTSSetJacobian()
512: @*/
513: PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
514: {
515:   DMTS           tsdm;

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

524: /*@C
525:    DMTSSetTransientVariable - sets function to transform from state to transient variables

527:    Logically Collective

529:    Input Parameters:
530: +  dm - DM to be used with TS
531: .  tvar - a function that transforms to transient variables
532: -  ctx - a context for tvar

534:     Calling sequence of tvar:
535: $     PetscErrorCode tvar(TS ts,Vec p,Vec c,void *ctx);

537: +   ts - timestep context
538: .   p - input vector (primitive form)
539: .   c - output vector, transient variables (conservative form)
540: -   ctx - [optional] user-defined function context

542:    Level: advanced

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

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

553: .seealso: TSSetTransientVariable(), DMTSGetTransientVariable(), DMTSSetIFunction(), DMTSSetIJacobian()
554: @*/
555: PetscErrorCode DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void *ctx)
556: {
557:   DMTS           dmts;

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

566: /*@C
567:    DMTSGetTransientVariable - gets function to transform from state to transient variables

569:    Logically Collective

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

574:    Output Parameters:
575: +  tvar - a function that transforms to transient variables
576: -  ctx - a context for tvar

578:    Level: advanced

580: .seealso: DMTSSetTransientVariable(), DMTSGetIFunction(), DMTSGetIJacobian()
581: @*/
582: PetscErrorCode DMTSGetTransientVariable(DM dm,TSTransientVariable *tvar,void *ctx)
583: {
584:   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 Parameter:
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: {
611:   DMTS           tsdm;

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

620: /*@C
621:    DMTSSetSolutionFunction - set TS solution evaluation function

623:    Not Collective

625:    Input Parameters:
626: +  dm - DM to be used with TS
627: .  func - solution function evaluation routine
628: -  ctx - context for solution evaluation

630:     Calling sequence of f:
631: $     PetscErrorCode f(TS ts,PetscReal t,Vec u,void *ctx);

633: +   ts - timestep context
634: .   t - current timestep
635: .   u - output vector
636: -   ctx - [optional] user-defined function context

638:    Level: advanced

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

645: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSGetSolutionFunction()
646: @*/
647: PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
648: {
649:   DMTS           tsdm;

652:   DMGetDMTSWrite(dm,&tsdm);
653:   if (func) tsdm->ops->solution = func;
654:   if (ctx)  tsdm->solutionctx   = ctx;
655:   return 0;
656: }

658: /*@C
659:    DMTSSetForcingFunction - set TS forcing function evaluation function

661:    Not Collective

663:    Input Parameters:
664: +  dm - DM to be used with TS
665: .  f - forcing function evaluation routine
666: -  ctx - context for solution evaluation

668:     Calling sequence of func:
669: $     PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx);

671: +   ts - timestep context
672: .   t - current timestep
673: .   f - output vector
674: -   ctx - [optional] user-defined function context

676:    Level: advanced

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

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

690:   DMGetDMTSWrite(dm,&tsdm);
691:   if (f)    tsdm->ops->forcing = f;
692:   if (ctx)  tsdm->forcingctx   = ctx;
693:   return 0;
694: }

696: /*@C
697:    DMTSGetForcingFunction - get TS forcing function evaluation function

699:    Not Collective

701:    Input Parameter:
702: .   dm - DM to be used with TS

704:    Output Parameters:
705: +  f - forcing function evaluation function; see TSForcingFunction for details
706: -  ctx - context for solution evaluation

708:    Level: advanced

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

715: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
716: @*/
717: PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
718: {
719:   DMTS           tsdm;

722:   DMGetDMTSWrite(dm,&tsdm);
723:   if (f)   *f   = tsdm->ops->forcing;
724:   if (ctx) *ctx = tsdm->forcingctx;
725:   return 0;
726: }

728: /*@C
729:    DMTSGetRHSFunction - get TS explicit residual evaluation function

731:    Not Collective

733:    Input Parameter:
734: .  dm - DM to be used with TS

736:    Output Parameters:
737: +  func - residual evaluation function, see TSSetRHSFunction() for calling sequence
738: -  ctx - context for residual evaluation

740:    Level: advanced

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

746: .seealso: DMTSSetContext(), DMTSSetRHSFunction(), TSSetRHSFunction()
747: @*/
748: PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
749: {
750:   DMTS           tsdm;

753:   DMGetDMTS(dm,&tsdm);
754:   if (func) *func = tsdm->ops->rhsfunction;
755:   if (ctx)  *ctx = tsdm->rhsfunctionctx;
756:   return 0;
757: }

759: /*@C
760:    DMTSSetIJacobian - set TS Jacobian evaluation function

762:    Not Collective

764:    Input Parameters:
765: +  dm - DM to be used with TS
766: .  func - Jacobian evaluation routine
767: -  ctx - context for residual evaluation

769:    Calling sequence of f:
770: $    PetscErrorCode f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx);

772: +  t    - time at step/stage being solved
773: .  U    - state vector
774: .  U_t  - time derivative of state vector
775: .  a    - shift
776: .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
777: .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
778: -  ctx  - [optional] user-defined context for matrix evaluation routine

780:    Level: advanced

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

787: .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSGetJacobian(), TSSetIJacobian(), TSSetIFunction()
788: @*/
789: PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
790: {
791:   DMTS           sdm;

794:   DMGetDMTSWrite(dm,&sdm);
795:   if (func) sdm->ops->ijacobian = func;
796:   if (ctx)  sdm->ijacobianctx   = ctx;
797:   return 0;
798: }

800: /*@C
801:    DMTSGetIJacobian - get TS Jacobian evaluation function

803:    Not Collective

805:    Input Parameter:
806: .  dm - DM to be used with TS

808:    Output Parameters:
809: +  func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
810: -  ctx - context for residual evaluation

812:    Level: advanced

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

819: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
820: @*/
821: PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
822: {
823:   DMTS           tsdm;

826:   DMGetDMTS(dm,&tsdm);
827:   if (func) *func = tsdm->ops->ijacobian;
828:   if (ctx)  *ctx = tsdm->ijacobianctx;
829:   return 0;
830: }

832: /*@C
833:    DMTSSetRHSJacobian - set TS Jacobian evaluation function

835:    Not Collective

837:    Input Parameters:
838: +  dm - DM to be used with TS
839: .  func - Jacobian evaluation routine
840: -  ctx - context for residual evaluation

842:    Calling sequence of func:
843: $     PetscErrorCode func(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);

845: +  t - current timestep
846: .  u - input vector
847: .  Amat - (approximate) Jacobian matrix
848: .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
849: -  ctx - [optional] user-defined context for matrix evaluation routine

851:    Level: advanced

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

858: .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetRHSJacobian()
859: @*/
860: PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
861: {
862:   DMTS           tsdm;

865:   DMGetDMTSWrite(dm,&tsdm);
866:   if (func) tsdm->ops->rhsjacobian = func;
867:   if (ctx)  tsdm->rhsjacobianctx = ctx;
868:   return 0;
869: }

871: /*@C
872:    DMTSGetRHSJacobian - get TS Jacobian evaluation function

874:    Not Collective

876:    Input Parameter:
877: .  dm - DM to be used with TS

879:    Output Parameters:
880: +  func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
881: -  ctx - context for residual evaluation

883:    Level: advanced

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

890: .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSSetRHSJacobian(), TSSetRHSJacobian()
891: @*/
892: PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
893: {
894:   DMTS           tsdm;

897:   DMGetDMTS(dm,&tsdm);
898:   if (func) *func = tsdm->ops->rhsjacobian;
899:   if (ctx)  *ctx = tsdm->rhsjacobianctx;
900:   return 0;
901: }

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

906:    Not Collective

908:    Input Parameters:
909: +  dm - DM to be used with TS
910: .  view - viewer function
911: -  load - loading function

913:    Level: advanced

915: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
916: @*/
917: PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
918: {
919:   DMTS           tsdm;

922:   DMGetDMTSWrite(dm,&tsdm);
923:   tsdm->ops->ifunctionview = view;
924:   tsdm->ops->ifunctionload = load;
925:   return 0;
926: }

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

931:    Not Collective

933:    Input Parameters:
934: +  dm - DM to be used with TS
935: .  view - viewer function
936: -  load - loading function

938:    Level: advanced

940: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
941: @*/
942: PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
943: {
944:   DMTS           tsdm;

947:   DMGetDMTSWrite(dm,&tsdm);
948:   tsdm->ops->ijacobianview = view;
949:   tsdm->ops->ijacobianload = load;
950:   return 0;
951: }