Actual source code: theta.c
petsc-3.3-p1 2012-06-15
1: /*
2: Code for timestepping with implicit Theta method
3: */
4: #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/
5: #include <petscsnesfas.h>
7: typedef struct {
8: Vec X,Xdot; /* Storage for one stage */
9: Vec affine; /* Affine vector needed for residual at beginning of step */
10: PetscBool extrapolate;
11: PetscBool endpoint;
12: PetscReal Theta;
13: PetscReal shift;
14: PetscReal stage_time;
15: } TS_Theta;
19: static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
20: {
21: TS_Theta *th = (TS_Theta*)ts->data;
25: if (X0) {
26: if (dm && dm != ts->dm) {
27: PetscObjectQuery((PetscObject)dm,"TSTheta_X0",(PetscObject*)X0);
28: if (!*X0) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"TSTheta_X0 has not been composed with DM from SNES");
29: } else *X0 = ts->vec_sol;
30: }
31: if (Xdot) {
32: if (dm && dm != ts->dm) {
33: PetscObjectQuery((PetscObject)dm,"TSTheta_Xdot",(PetscObject*)Xdot);
34: if (!*Xdot) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"TSTheta_Xdot has not been composed with DM from SNES");
35: } else *Xdot = th->Xdot;
36: }
37: return(0);
38: }
42: static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
43: {
44: Vec X0,Xdot;
48: DMCreateGlobalVector(coarse,&X0);
49: DMCreateGlobalVector(coarse,&Xdot);
50: /* Oh noes, this would create a loop because the Vec holds a reference to the DM.
51: Making a PetscContainer to hold these Vecs would make the following call succeed, but would create a reference loop.
52: Need to decide on a way to break the reference counting loop.
53: */
54: PetscObjectCompose((PetscObject)coarse,"TSTheta_X0",(PetscObject)X0);
55: PetscObjectCompose((PetscObject)coarse,"TSTheta_Xdot",(PetscObject)Xdot);
56: VecDestroy(&X0);
57: VecDestroy(&Xdot);
58: return(0);
59: }
63: static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
64: {
65: TS ts = (TS)ctx;
67: Vec X0,Xdot,X0_c,Xdot_c;
70: TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);
71: TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);
72: MatRestrict(restrct,X0,X0_c);
73: MatRestrict(restrct,Xdot,Xdot_c);
74: VecPointwiseMult(X0_c,rscale,X0_c);
75: VecPointwiseMult(Xdot_c,rscale,Xdot_c);
76: return(0);
77: }
81: static PetscErrorCode TSStep_Theta(TS ts)
82: {
83: TS_Theta *th = (TS_Theta*)ts->data;
84: PetscInt its,lits;
85: PetscReal next_time_step;
86: SNESConvergedReason snesreason;
87: PetscErrorCode ierr;
90: next_time_step = ts->time_step;
91: th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
92: th->shift = 1./(th->Theta*ts->time_step);
94: if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
95: VecZeroEntries(th->Xdot);
96: if (!th->affine) {VecDuplicate(ts->vec_sol,&th->affine);}
97: TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);
98: VecScale(th->affine,(th->Theta-1.)/th->Theta);
99: }
100: if (th->extrapolate) {
101: VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);
102: } else {
103: VecCopy(ts->vec_sol,th->X);
104: }
105: SNESSolve(ts->snes,th->affine,th->X);
106: SNESGetIterationNumber(ts->snes,&its);
107: SNESGetLinearSolveIterations(ts->snes,&lits);
108: SNESGetConvergedReason(ts->snes,&snesreason);
109: ts->snes_its += its; ts->ksp_its += lits;
110: if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
111: ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
112: PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
113: return(0);
114: }
115: if (th->endpoint) {
116: VecCopy(th->X,ts->vec_sol);
117: } else {
118: VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);
119: VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);
120: }
121: ts->ptime += ts->time_step;
122: ts->time_step = next_time_step;
123: ts->steps++;
124: return(0);
125: }
129: static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
130: {
131: TS_Theta *th = (TS_Theta*)ts->data;
132: PetscReal alpha = t - ts->ptime;
136: VecCopy(ts->vec_sol,th->X);
137: if (th->endpoint) alpha *= th->Theta;
138: VecWAXPY(X,alpha,th->Xdot,th->X);
139: return(0);
140: }
142: /*------------------------------------------------------------*/
145: static PetscErrorCode TSReset_Theta(TS ts)
146: {
147: TS_Theta *th = (TS_Theta*)ts->data;
148: PetscErrorCode ierr;
151: VecDestroy(&th->X);
152: VecDestroy(&th->Xdot);
153: VecDestroy(&th->affine);
154: return(0);
155: }
159: static PetscErrorCode TSDestroy_Theta(TS ts)
160: {
161: PetscErrorCode ierr;
164: TSReset_Theta(ts);
165: PetscFree(ts->data);
166: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);
167: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);
168: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);
169: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);
170: return(0);
171: }
173: /*
174: This defines the nonlinear equation that is to be solved with SNES
175: G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
176: */
179: static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
180: {
181: TS_Theta *th = (TS_Theta*)ts->data;
183: Vec X0,Xdot;
184: DM dm,dmsave;
187: SNESGetDM(snes,&dm);
188: /* When using the endpoint variant, this is actually 1/Theta * Xdot */
189: TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);
190: VecAXPBYPCZ(Xdot,-th->shift,th->shift,0,X0,x);
192: /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
193: dmsave = ts->dm;
194: ts->dm = dm;
195: TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);
196: ts->dm = dmsave;
197: return(0);
198: }
202: static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
203: {
204: TS_Theta *th = (TS_Theta*)ts->data;
206: Vec Xdot;
207: DM dm,dmsave;
210: SNESGetDM(snes,&dm);
212: /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
213: TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);
215: dmsave = ts->dm;
216: ts->dm = dm;
217: TSComputeIJacobian(ts,th->stage_time,x,Xdot,th->shift,A,B,str,PETSC_FALSE);
218: ts->dm = dmsave;
219: return(0);
220: }
224: static PetscErrorCode TSSetUp_Theta(TS ts)
225: {
226: TS_Theta *th = (TS_Theta*)ts->data;
228: SNES snes;
229: DM dm;
232: VecDuplicate(ts->vec_sol,&th->X);
233: VecDuplicate(ts->vec_sol,&th->Xdot);
234: TSGetSNES(ts,&snes);
235: TSGetDM(ts,&dm);
236: if (dm) {
237: DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);
238: }
239: return(0);
240: }
241: /*------------------------------------------------------------*/
245: static PetscErrorCode TSSetFromOptions_Theta(TS ts)
246: {
247: TS_Theta *th = (TS_Theta*)ts->data;
251: PetscOptionsHead("Theta ODE solver options");
252: {
253: PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);
254: PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);
255: PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,PETSC_NULL);
256: SNESSetFromOptions(ts->snes);
257: }
258: PetscOptionsTail();
259: return(0);
260: }
264: static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
265: {
266: TS_Theta *th = (TS_Theta*)ts->data;
267: PetscBool iascii;
268: PetscErrorCode ierr;
271: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
272: if (iascii) {
273: PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);
274: PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");
275: }
276: SNESView(ts->snes,viewer);
277: return(0);
278: }
280: EXTERN_C_BEGIN
283: PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
284: {
285: TS_Theta *th = (TS_Theta*)ts->data;
288: *theta = th->Theta;
289: return(0);
290: }
294: PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
295: {
296: TS_Theta *th = (TS_Theta*)ts->data;
299: if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
300: th->Theta = theta;
301: return(0);
302: }
306: PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
307: {
308: TS_Theta *th = (TS_Theta*)ts->data;
311: *endpoint = th->endpoint;
312: return(0);
313: }
317: PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
318: {
319: TS_Theta *th = (TS_Theta*)ts->data;
322: th->endpoint = flg;
323: return(0);
324: }
325: EXTERN_C_END
327: /* ------------------------------------------------------------ */
328: /*MC
329: TSTHETA - DAE solver using the implicit Theta method
331: Level: beginner
333: Notes:
334: This method can be applied to DAE.
336: This method is cast as a 1-stage implicit Runge-Kutta method.
338: .vb
339: Theta | Theta
340: -------------
341: | 1
342: .ve
344: For the default Theta=0.5, this is also known as the implicit midpoint rule.
346: When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
348: .vb
349: 0 | 0 0
350: 1 | 1-Theta Theta
351: -------------------
352: | 1-Theta Theta
353: .ve
355: For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
357: To apply a diagonally implicit RK method to DAE, the stage formula
359: $ Y_i = X + h sum_j a_ij Y'_j
361: is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i)
363: .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
365: M*/
366: EXTERN_C_BEGIN
369: PetscErrorCode TSCreate_Theta(TS ts)
370: {
371: TS_Theta *th;
375: ts->ops->reset = TSReset_Theta;
376: ts->ops->destroy = TSDestroy_Theta;
377: ts->ops->view = TSView_Theta;
378: ts->ops->setup = TSSetUp_Theta;
379: ts->ops->step = TSStep_Theta;
380: ts->ops->interpolate = TSInterpolate_Theta;
381: ts->ops->setfromoptions = TSSetFromOptions_Theta;
382: ts->ops->snesfunction = SNESTSFormFunction_Theta;
383: ts->ops->snesjacobian = SNESTSFormJacobian_Theta;
385: PetscNewLog(ts,TS_Theta,&th);
386: ts->data = (void*)th;
388: th->extrapolate = PETSC_FALSE;
389: th->Theta = 0.5;
391: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);
392: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);
393: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);
394: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);
395: return(0);
396: }
397: EXTERN_C_END
401: /*@
402: TSThetaGetTheta - Get the abscissa of the stage in (0,1].
404: Not Collective
406: Input Parameter:
407: . ts - timestepping context
409: Output Parameter:
410: . theta - stage abscissa
412: Note:
413: Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
415: Level: Advanced
417: .seealso: TSThetaSetTheta()
418: @*/
419: PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta)
420: {
426: PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));
427: return(0);
428: }
432: /*@
433: TSThetaSetTheta - Set the abscissa of the stage in (0,1].
435: Not Collective
437: Input Parameter:
438: + ts - timestepping context
439: - theta - stage abscissa
441: Options Database:
442: . -ts_theta_theta <theta>
444: Level: Intermediate
446: .seealso: TSThetaGetTheta()
447: @*/
448: PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta)
449: {
454: PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));
455: return(0);
456: }
460: /*@
461: TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
463: Not Collective
465: Input Parameter:
466: . ts - timestepping context
468: Output Parameter:
469: . endpoint - PETSC_TRUE when using the endpoint variant
471: Level: Advanced
473: .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
474: @*/
475: PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
476: {
482: PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));
483: return(0);
484: }
488: /*@
489: TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
491: Not Collective
493: Input Parameter:
494: + ts - timestepping context
495: - flg - PETSC_TRUE to use the endpoint variant
497: Options Database:
498: . -ts_theta_endpoint <flg>
500: Level: Intermediate
502: .seealso: TSTHETA, TSCN
503: @*/
504: PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
505: {
510: PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));
511: return(0);
512: }
514: /*
515: * TSBEULER and TSCN are straightforward specializations of TSTHETA.
516: * The creation functions for these specializations are below.
517: */
521: static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
522: {
526: SNESView(ts->snes,viewer);
527: return(0);
528: }
530: /*MC
531: TSBEULER - ODE solver using the implicit backward Euler method
533: Level: beginner
535: .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
537: M*/
538: EXTERN_C_BEGIN
541: PetscErrorCode TSCreate_BEuler(TS ts)
542: {
546: TSCreate_Theta(ts);
547: TSThetaSetTheta(ts,1.0);
548: ts->ops->view = TSView_BEuler;
549: return(0);
550: }
551: EXTERN_C_END
555: static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
556: {
560: SNESView(ts->snes,viewer);
561: return(0);
562: }
564: /*MC
565: TSCN - ODE solver using the implicit Crank-Nicolson method.
567: Level: beginner
569: Notes:
570: TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
572: $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
574: .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
576: M*/
577: EXTERN_C_BEGIN
580: PetscErrorCode TSCreate_CN(TS ts)
581: {
585: TSCreate_Theta(ts);
586: TSThetaSetTheta(ts,0.5);
587: TSThetaSetEndpoint(ts,PETSC_TRUE);
588: ts->ops->view = TSView_CN;
589: return(0);
590: }
591: EXTERN_C_END