Actual source code: dmts.c
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 = NULL; return(0);}
12: if ((*kdm)->ops->destroy) {((*kdm)->ops->destroy)(*kdm);}
13: PetscHeaderDestroy(kdm);
14: return(0);
15: }
17: PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer)
18: {
22: PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,NULL,PETSC_FUNCTION);
23: PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,NULL,PETSC_FUNCTION);
24: PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,NULL,PETSC_FUNCTION);
25: if (kdm->ops->ifunctionload) {
26: (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);
27: }
28: PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,NULL,PETSC_FUNCTION);
29: PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,NULL,PETSC_FUNCTION);
30: PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,NULL,PETSC_FUNCTION);
31: if (kdm->ops->ijacobianload) {
32: (*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer);
33: }
34: return(0);
35: }
37: PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer)
38: {
40: PetscBool isascii,isbinary;
43: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
44: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
45: if (isascii) {
46: #if defined(PETSC_SERIALIZE_FUNCTIONS)
47: const char *fname;
49: PetscFPTFind(kdm->ops->ifunction,&fname);
50: if (fname) {
51: PetscViewerASCIIPrintf(viewer," IFunction used by TS: %s\n",fname);
52: }
53: PetscFPTFind(kdm->ops->ijacobian,&fname);
54: if (fname) {
55: PetscViewerASCIIPrintf(viewer," IJacobian function used by TS: %s\n",fname);
56: }
57: #endif
58: } else if (isbinary) {
59: struct {
60: TSIFunction ifunction;
61: } funcstruct;
62: struct {
63: PetscErrorCode (*ifunctionview)(void*,PetscViewer);
64: } funcviewstruct;
65: struct {
66: PetscErrorCode (*ifunctionload)(void**,PetscViewer);
67: } funcloadstruct;
68: struct {
69: TSIJacobian ijacobian;
70: } jacstruct;
71: struct {
72: PetscErrorCode (*ijacobianview)(void*,PetscViewer);
73: } jacviewstruct;
74: struct {
75: PetscErrorCode (*ijacobianload)(void**,PetscViewer);
76: } jacloadstruct;
78: funcstruct.ifunction = kdm->ops->ifunction;
79: funcviewstruct.ifunctionview = kdm->ops->ifunctionview;
80: funcloadstruct.ifunctionload = kdm->ops->ifunctionload;
81: PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION);
82: PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION);
83: PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION);
84: if (kdm->ops->ifunctionview) {
85: (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);
86: }
87: jacstruct.ijacobian = kdm->ops->ijacobian;
88: jacviewstruct.ijacobianview = kdm->ops->ijacobianview;
89: jacloadstruct.ijacobianload = kdm->ops->ijacobianload;
90: PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION);
91: PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION);
92: PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION);
93: if (kdm->ops->ijacobianview) {
94: (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);
95: }
96: }
97: return(0);
98: }
100: static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm)
101: {
105: TSInitializePackage();
106: PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);
107: return(0);
108: }
110: /* Attaches the DMTS to the coarse level.
111: * Under what conditions should we copy versus duplicate?
112: */
113: static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx)
114: {
118: DMCopyDMTS(dm,dmc);
119: return(0);
120: }
122: /* This could restrict auxiliary information to the coarse level.
123: */
124: static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
125: {
128: return(0);
129: }
131: static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx)
132: {
136: DMCopyDMTS(dm,subdm);
137: return(0);
138: }
140: /* This could restrict auxiliary information to the coarse level.
141: */
142: static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
143: {
145: return(0);
146: }
148: /*@C
149: DMTSCopy - copies the information in a DMTS to another DMTS
151: Not Collective
153: Input Argument:
154: + kdm - Original DMTS
155: - nkdm - DMTS to receive the data, should have been created with DMTSCreate()
157: Level: developer
159: .seealso: DMTSCreate(), DMTSDestroy()
160: @*/
161: PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm)
162: {
168: nkdm->ops->rhsfunction = kdm->ops->rhsfunction;
169: nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian;
170: nkdm->ops->ifunction = kdm->ops->ifunction;
171: nkdm->ops->ijacobian = kdm->ops->ijacobian;
172: nkdm->ops->i2function = kdm->ops->i2function;
173: nkdm->ops->i2jacobian = kdm->ops->i2jacobian;
174: nkdm->ops->solution = kdm->ops->solution;
175: nkdm->ops->destroy = kdm->ops->destroy;
176: nkdm->ops->duplicate = kdm->ops->duplicate;
178: nkdm->rhsfunctionctx = kdm->rhsfunctionctx;
179: nkdm->rhsjacobianctx = kdm->rhsjacobianctx;
180: nkdm->ifunctionctx = kdm->ifunctionctx;
181: nkdm->ijacobianctx = kdm->ijacobianctx;
182: nkdm->i2functionctx = kdm->i2functionctx;
183: nkdm->i2jacobianctx = kdm->i2jacobianctx;
184: nkdm->solutionctx = kdm->solutionctx;
186: nkdm->data = kdm->data;
188: /*
189: nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
190: nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
191: nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
192: */
194: /* implementation specific copy hooks */
195: if (kdm->ops->duplicate) {(*kdm->ops->duplicate)(kdm,nkdm);}
196: return(0);
197: }
199: /*@C
200: DMGetDMTS - get read-only private DMTS context from a DM
202: Not Collective
204: Input Argument:
205: . dm - DM to be used with TS
207: Output Argument:
208: . tsdm - private DMTS context
210: Level: developer
212: Notes:
213: Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible.
215: .seealso: DMGetDMTSWrite()
216: @*/
217: PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm)
218: {
223: *tsdm = (DMTS) dm->dmts;
224: if (!*tsdm) {
225: PetscInfo(dm,"Creating new DMTS\n");
226: DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm);
227: dm->dmts = (PetscObject) *tsdm;
228: (*tsdm)->originaldm = dm;
229: DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);
230: DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);
231: }
232: return(0);
233: }
235: /*@C
236: DMGetDMTSWrite - get write access to private DMTS context from a DM
238: Not Collective
240: Input Argument:
241: . dm - DM to be used with TS
243: Output Argument:
244: . tsdm - private DMTS context
246: Level: developer
248: .seealso: DMGetDMTS()
249: @*/
250: PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm)
251: {
253: DMTS sdm;
257: DMGetDMTS(dm,&sdm);
258: if (!sdm->originaldm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DMTS has a NULL originaldm");
259: if (sdm->originaldm != dm) { /* Copy on write */
260: DMTS oldsdm = sdm;
261: PetscInfo(dm,"Copying DMTS due to write\n");
262: DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);
263: DMTSCopy(oldsdm,sdm);
264: DMTSDestroy((DMTS*)&dm->dmts);
265: dm->dmts = (PetscObject) sdm;
266: sdm->originaldm = dm;
267: }
268: *tsdm = sdm;
269: return(0);
270: }
272: /*@C
273: DMCopyDMTS - copies a DM context to a new DM
275: Logically Collective
277: Input Arguments:
278: + dmsrc - DM to obtain context from
279: - dmdest - DM to add context to
281: Level: developer
283: Note:
284: The context is copied by reference. This function does not ensure that a context exists.
286: .seealso: DMGetDMTS(), TSSetDM()
287: @*/
288: PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest)
289: {
295: DMTSDestroy((DMTS*)&dmdest->dmts);
296: dmdest->dmts = dmsrc->dmts;
297: PetscObjectReference(dmdest->dmts);
298: DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);
299: DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);
300: return(0);
301: }
303: /*@C
304: DMTSSetIFunction - set TS implicit function evaluation function
306: Not Collective
308: Input Arguments:
309: + dm - DM to be used with TS
310: . func - function evaluating f(t,u,u_t)
311: - ctx - context for residual evaluation
313: Calling sequence of func:
314: $ PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
316: + t - time at step/stage being solved
317: . u - state vector
318: . u_t - time derivative of state vector
319: . F - function vector
320: - ctx - [optional] user-defined context for matrix evaluation routine
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(), TSSetIFunction(), 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: }
344: /*@C
345: DMTSGetIFunction - get TS implicit residual evaluation function
347: Not Collective
349: Input Argument:
350: . dm - DM to be used with TS
352: Output Arguments:
353: + func - function evaluation function, see TSSetIFunction() for calling sequence
354: - ctx - context for residual evaluation
356: Level: advanced
358: Note:
359: TSGetFunction() is normally used, but it calls this function internally because the user context is actually
360: associated with the DM.
362: .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction()
363: @*/
364: PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx)
365: {
367: DMTS tsdm;
371: DMGetDMTS(dm,&tsdm);
372: if (func) *func = tsdm->ops->ifunction;
373: if (ctx) *ctx = tsdm->ifunctionctx;
374: return(0);
375: }
377: /*@C
378: DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems
380: Not Collective
382: Input Arguments:
383: + dm - DM to be used with TS
384: . fun - function evaluation routine
385: - ctx - context for residual evaluation
387: Calling sequence of fun:
388: $ PetscErrorCode fun(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,Vec F,ctx);
390: + t - time at step/stage being solved
391: . U - state vector
392: . U_t - time derivative of state vector
393: . U_tt - second time derivative of state vector
394: . F - function vector
395: - ctx - [optional] user-defined context for matrix evaluation routine (may be NULL)
397: Level: advanced
399: Note:
400: TSSetI2Function() is normally used, but it calls this function internally because the user context is actually
401: associated with the DM.
403: .seealso: TSSetI2Function()
404: @*/
405: PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx)
406: {
407: DMTS tsdm;
412: DMGetDMTSWrite(dm,&tsdm);
413: if (fun) tsdm->ops->i2function = fun;
414: if (ctx) tsdm->i2functionctx = ctx;
415: return(0);
416: }
418: /*@C
419: DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems
421: Not Collective
423: Input Argument:
424: . dm - DM to be used with TS
426: Output Arguments:
427: + fun - function evaluation function, see TSSetI2Function() for calling sequence
428: - ctx - context for residual evaluation
430: Level: advanced
432: Note:
433: TSGetI2Function() is normally used, but it calls this function internally because the user context is actually
434: associated with the DM.
436: .seealso: DMTSSetI2Function(),TSGetI2Function()
437: @*/
438: PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx)
439: {
440: DMTS tsdm;
445: DMGetDMTS(dm,&tsdm);
446: if (fun) *fun = tsdm->ops->i2function;
447: if (ctx) *ctx = tsdm->i2functionctx;
448: return(0);
449: }
451: /*@C
452: DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems
454: Not Collective
456: Input Arguments:
457: + dm - DM to be used with TS
458: . fun - Jacobian evaluation routine
459: - ctx - context for Jacobian evaluation
461: Calling sequence of jac:
462: $ PetscErrorCode jac(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,PetscReal v,PetscReal a,Mat J,Mat P,void *ctx);
464: + t - time at step/stage being solved
465: . U - state vector
466: . U_t - time derivative of state vector
467: . U_tt - second time derivative of state vector
468: . v - shift for U_t
469: . a - shift for U_tt
470: . 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
471: . P - preconditioning matrix for J, may be same as J
472: - ctx - [optional] user-defined context for matrix evaluation routine
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: }
495: /*@C
496: DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems
498: Not Collective
500: Input Argument:
501: . dm - DM to be used with TS
503: Output Arguments:
504: + jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence
505: - ctx - context for Jacobian evaluation
507: Level: advanced
509: Note:
510: TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually
511: associated with the DM.
513: .seealso: DMTSSetI2Jacobian(),TSGetI2Jacobian()
514: @*/
515: PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx)
516: {
517: DMTS tsdm;
522: DMGetDMTS(dm,&tsdm);
523: if (jac) *jac = tsdm->ops->i2jacobian;
524: if (ctx) *ctx = tsdm->i2jacobianctx;
525: return(0);
526: }
528: /*@C
529: DMTSSetRHSFunction - set TS explicit residual evaluation function
531: Not Collective
533: Input Arguments:
534: + dm - DM to be used with TS
535: . func - RHS function evaluation routine
536: - ctx - context for residual evaluation
538: Calling sequence of func:
539: $ PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec F,void *ctx);
541: + ts - timestep context
542: . t - current timestep
543: . u - input vector
544: . F - function vector
545: - ctx - [optional] user-defined function context
547: Level: advanced
549: Note:
550: TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually
551: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
552: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
554: .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSSetJacobian()
555: @*/
556: PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx)
557: {
559: DMTS tsdm;
563: DMGetDMTSWrite(dm,&tsdm);
564: if (func) tsdm->ops->rhsfunction = func;
565: if (ctx) tsdm->rhsfunctionctx = ctx;
566: return(0);
567: }
569: /*@C
570: DMTSSetTransientVariable - sets function to transform from state to transient variables
572: Logically Collective
574: Input Arguments:
575: + dm - DM to be used with TS
576: . tvar - a function that transforms to transient variables
577: - ctx - a context for tvar
579: Calling sequence of tvar:
580: $ PetscErrorCode tvar(TS ts,Vec p,Vec c,void *ctx);
582: + ts - timestep context
583: . p - input vector (primative form)
584: . c - output vector, transient variables (conservative form)
585: - ctx - [optional] user-defined function context
587: Level: advanced
589: Notes:
590: This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF)
591: can be conservative. In this context, primitive variables P are used to model the state (e.g., because they lead to
592: well-conditioned formulations even in limiting cases such as low-Mach or zero porosity). The transient variable is
593: C(P), specified by calling this function. An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
594: evaluated via the chain rule, as in
596: dF/dP + shift * dF/dCdot dC/dP.
598: .seealso: TSSetTransientVariable(), DMTSGetTransientVariable(), DMTSSetIFunction(), DMTSSetIJacobian()
599: @*/
600: PetscErrorCode DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void *ctx)
601: {
603: DMTS dmts;
607: DMGetDMTSWrite(dm,&dmts);
608: dmts->ops->transientvar = tvar;
609: dmts->transientvarctx = ctx;
610: return(0);
611: }
613: /*@C
614: DMTSGetTransientVariable - gets function to transform from state to transient variables
616: Logically Collective
618: Input Arguments:
619: . dm - DM to be used with TS
621: Output Arguments:
622: + tvar - a function that transforms to transient variables
623: - ctx - a context for tvar
625: Level: advanced
627: .seealso: DMTSSetTransientVariable(), DMTSGetIFunction(), DMTSGetIJacobian()
628: @*/
629: PetscErrorCode DMTSGetTransientVariable(DM dm,TSTransientVariable *tvar,void *ctx)
630: {
632: DMTS dmts;
636: DMGetDMTS(dm,&dmts);
637: if (tvar) *tvar = dmts->ops->transientvar;
638: if (ctx) *(void**)ctx = dmts->transientvarctx;
639: return(0);
640: }
642: /*@C
643: DMTSGetSolutionFunction - gets the TS solution evaluation function
645: Not Collective
647: Input Arguments:
648: . dm - DM to be used with TS
650: Output Parameters:
651: + func - solution function evaluation function, see TSSetSolution() for calling sequence
652: - ctx - context for solution evaluation
654: Level: advanced
656: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSSetSolutionFunction()
657: @*/
658: PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx)
659: {
661: DMTS tsdm;
665: DMGetDMTS(dm,&tsdm);
666: if (func) *func = tsdm->ops->solution;
667: if (ctx) *ctx = tsdm->solutionctx;
668: return(0);
669: }
671: /*@C
672: DMTSSetSolutionFunction - set TS solution evaluation function
674: Not Collective
676: Input Arguments:
677: + dm - DM to be used with TS
678: . func - solution function evaluation routine
679: - ctx - context for solution evaluation
681: Calling sequence of f:
682: $ PetscErrorCode f(TS ts,PetscReal t,Vec u,void *ctx);
684: + ts - timestep context
685: . t - current timestep
686: . u - output vector
687: - ctx - [optional] user-defined function context
689: Level: advanced
691: Note:
692: TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually
693: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
694: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
696: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSGetSolutionFunction()
697: @*/
698: PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx)
699: {
701: DMTS tsdm;
705: DMGetDMTSWrite(dm,&tsdm);
706: if (func) tsdm->ops->solution = func;
707: if (ctx) tsdm->solutionctx = ctx;
708: return(0);
709: }
711: /*@C
712: DMTSSetForcingFunction - set TS forcing function evaluation function
714: Not Collective
716: Input Arguments:
717: + dm - DM to be used with TS
718: . f - forcing function evaluation routine
719: - ctx - context for solution evaluation
721: Calling sequence of func:
722: $ PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx);
724: + ts - timestep context
725: . t - current timestep
726: . f - output vector
727: - ctx - [optional] user-defined function context
729: Level: advanced
731: Note:
732: TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually
733: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
734: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
736: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
737: @*/
738: PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx)
739: {
741: DMTS tsdm;
745: DMGetDMTSWrite(dm,&tsdm);
746: if (f) tsdm->ops->forcing = f;
747: if (ctx) tsdm->forcingctx = ctx;
748: return(0);
749: }
752: /*@C
753: DMTSGetForcingFunction - get TS forcing function evaluation function
755: Not Collective
757: Input Argument:
758: . dm - DM to be used with TS
760: Output Arguments:
761: + f - forcing function evaluation function; see TSForcingFunction for details
762: - ctx - context for solution evaluation
764: Level: advanced
766: Note:
767: TSSetForcingFunction() 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 residual.
771: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction()
772: @*/
773: PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx)
774: {
776: DMTS tsdm;
780: DMGetDMTSWrite(dm,&tsdm);
781: if (f) *f = tsdm->ops->forcing;
782: if (ctx) *ctx = tsdm->forcingctx;
783: return(0);
784: }
786: /*@C
787: DMTSGetRHSFunction - get TS explicit residual evaluation function
789: Not Collective
791: Input Argument:
792: . dm - DM to be used with TS
794: Output Arguments:
795: + func - residual evaluation function, see TSSetRHSFunction() for calling sequence
796: - ctx - context for residual evaluation
798: Level: advanced
800: Note:
801: TSGetFunction() is normally used, but it calls this function internally because the user context is actually
802: associated with the DM.
804: .seealso: DMTSSetContext(), DMTSSetRHSFunction(), TSSetRHSFunction()
805: @*/
806: PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx)
807: {
809: DMTS tsdm;
813: DMGetDMTS(dm,&tsdm);
814: if (func) *func = tsdm->ops->rhsfunction;
815: if (ctx) *ctx = tsdm->rhsfunctionctx;
816: return(0);
817: }
819: /*@C
820: DMTSSetIJacobian - set TS Jacobian evaluation function
822: Not Collective
824: Input Argument:
825: + dm - DM to be used with TS
826: . func - Jacobian evaluation routine
827: - ctx - context for residual evaluation
829: Calling sequence of f:
830: $ PetscErrorCode f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx);
832: + t - time at step/stage being solved
833: . U - state vector
834: . U_t - time derivative of state vector
835: . a - shift
836: . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
837: . Pmat - matrix used for constructing preconditioner, usually the same as Amat
838: - ctx - [optional] user-defined context for matrix evaluation routine
840: Level: advanced
842: Note:
843: TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
844: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
845: not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
847: .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSGetJacobian(), TSSetIJacobian(), TSSetIFunction()
848: @*/
849: PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx)
850: {
852: DMTS sdm;
856: DMGetDMTSWrite(dm,&sdm);
857: if (func) sdm->ops->ijacobian = func;
858: if (ctx) sdm->ijacobianctx = ctx;
859: return(0);
860: }
862: /*@C
863: DMTSGetIJacobian - get TS Jacobian evaluation function
865: Not Collective
867: Input Argument:
868: . dm - DM to be used with TS
870: Output Arguments:
871: + func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence
872: - ctx - context for residual evaluation
874: Level: advanced
876: Note:
877: TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
878: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
879: not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
881: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
882: @*/
883: PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx)
884: {
886: DMTS tsdm;
890: DMGetDMTS(dm,&tsdm);
891: if (func) *func = tsdm->ops->ijacobian;
892: if (ctx) *ctx = tsdm->ijacobianctx;
893: return(0);
894: }
896: /*@C
897: DMTSSetRHSJacobian - set TS Jacobian evaluation function
899: Not Collective
901: Input Argument:
902: + dm - DM to be used with TS
903: . func - Jacobian evaluation routine
904: - ctx - context for residual evaluation
906: Calling sequence of func:
907: $ PetscErrorCode func(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);
909: + t - current timestep
910: . u - input vector
911: . Amat - (approximate) Jacobian matrix
912: . Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
913: - ctx - [optional] user-defined context for matrix evaluation routine
915: Level: advanced
917: Note:
918: TSSetJacobian() is normally used, but it calls this function internally because the user context is actually
919: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
920: not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
922: .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetRHSJacobian()
923: @*/
924: PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx)
925: {
927: DMTS tsdm;
931: DMGetDMTSWrite(dm,&tsdm);
932: if (func) tsdm->ops->rhsjacobian = func;
933: if (ctx) tsdm->rhsjacobianctx = ctx;
934: return(0);
935: }
937: /*@C
938: DMTSGetRHSJacobian - get TS Jacobian evaluation function
940: Not Collective
942: Input Argument:
943: . dm - DM to be used with TS
945: Output Arguments:
946: + func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence
947: - ctx - context for residual evaluation
949: Level: advanced
951: Note:
952: TSGetJacobian() is normally used, but it calls this function internally because the user context is actually
953: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
954: not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
956: .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSSetRHSJacobian(), TSSetRHSJacobian()
957: @*/
958: PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx)
959: {
961: DMTS tsdm;
965: DMGetDMTS(dm,&tsdm);
966: if (func) *func = tsdm->ops->rhsjacobian;
967: if (ctx) *ctx = tsdm->rhsjacobianctx;
968: return(0);
969: }
971: /*@C
972: DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context
974: Not Collective
976: Input Arguments:
977: + dm - DM to be used with TS
978: . view - viewer function
979: - load - loading function
981: Level: advanced
983: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
984: @*/
985: PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
986: {
988: DMTS tsdm;
992: DMGetDMTSWrite(dm,&tsdm);
993: tsdm->ops->ifunctionview = view;
994: tsdm->ops->ifunctionload = load;
995: return(0);
996: }
998: /*@C
999: DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context
1001: Not Collective
1003: Input Arguments:
1004: + dm - DM to be used with TS
1005: . view - viewer function
1006: - load - loading function
1008: Level: advanced
1010: .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian()
1011: @*/
1012: PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer))
1013: {
1015: DMTS tsdm;
1019: DMGetDMTSWrite(dm,&tsdm);
1020: tsdm->ops->ijacobianview = view;
1021: tsdm->ops->ijacobianload = load;
1022: return(0);
1023: }