Actual source code: theta.c
petsc-3.6.0 2015-06-09
1: /*
2: Code for timestepping with implicit Theta method
3: */
4: #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
5: #include <petscsnes.h>
6: #include <petscdm.h>
7: #include <petscmat.h>
9: typedef struct {
10: Vec X,Xdot; /* Storage for one stage */
11: Vec X0; /* work vector to store X0 */
12: Vec affine; /* Affine vector needed for residual at beginning of step */
13: Vec *VecDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage*/
14: Vec *VecDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage*/
15: Vec *VecSensiTemp; /* Vector to be timed with Jacobian transpose*/
16: PetscBool extrapolate;
17: PetscBool endpoint;
18: PetscReal Theta;
19: PetscReal stage_time;
20: TSStepStatus status;
21: char *name;
22: PetscInt order;
23: PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */
24: PetscBool adapt; /* use time-step adaptivity ? */
25: } TS_Theta;
29: static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
30: {
31: TS_Theta *th = (TS_Theta*)ts->data;
35: if (X0) {
36: if (dm && dm != ts->dm) {
37: DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);
38: } else *X0 = ts->vec_sol;
39: }
40: if (Xdot) {
41: if (dm && dm != ts->dm) {
42: DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);
43: } else *Xdot = th->Xdot;
44: }
45: return(0);
46: }
51: static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
52: {
56: if (X0) {
57: if (dm && dm != ts->dm) {
58: DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);
59: }
60: }
61: if (Xdot) {
62: if (dm && dm != ts->dm) {
63: DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);
64: }
65: }
66: return(0);
67: }
71: static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
72: {
75: return(0);
76: }
80: static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
81: {
82: TS ts = (TS)ctx;
84: Vec X0,Xdot,X0_c,Xdot_c;
87: TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);
88: TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);
89: MatRestrict(restrct,X0,X0_c);
90: MatRestrict(restrct,Xdot,Xdot_c);
91: VecPointwiseMult(X0_c,rscale,X0_c);
92: VecPointwiseMult(Xdot_c,rscale,Xdot_c);
93: TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);
94: TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);
95: return(0);
96: }
100: static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
101: {
104: return(0);
105: }
109: static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
110: {
111: TS ts = (TS)ctx;
113: Vec X0,Xdot,X0_sub,Xdot_sub;
116: TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);
117: TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);
119: VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);
120: VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);
122: VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);
123: VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);
125: TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);
126: TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);
127: return(0);
128: }
132: static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done)
133: {
135: TS_Theta *th = (TS_Theta*)ts->data;
138: if (order == 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"No time-step adaptivity implemented for 1st order theta method; Run with -ts_adapt_type none");
139: if (order == th->order) {
140: if (th->endpoint) {
141: VecCopy(th->X,U);
142: } else {
143: PetscReal shift = 1./(th->Theta*ts->time_step);
144: VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);
145: VecAXPY(U,ts->time_step,th->Xdot);
146: }
147: } else if (order == th->order-1 && order) {
148: VecWAXPY(U,ts->time_step,th->Xdot,th->X0);
149: }
150: return(0);
151: }
155: static PetscErrorCode TSRollBack_Theta(TS ts)
156: {
157: TS_Theta *th = (TS_Theta*)ts->data;
161: VecCopy(th->X0,ts->vec_sol);
162: th->status = TS_STEP_INCOMPLETE;
163: return(0);
164: }
168: static PetscErrorCode TSStep_Theta(TS ts)
169: {
170: TS_Theta *th = (TS_Theta*)ts->data;
171: PetscInt its,lits,reject,next_scheme;
172: PetscReal next_time_step;
173: TSAdapt adapt;
174: PetscBool stageok,accept = PETSC_TRUE;
178: th->status = TS_STEP_INCOMPLETE;
179: VecCopy(ts->vec_sol,th->X0);
180: for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) {
181: PetscReal shift = 1./(th->Theta*ts->time_step);
182: th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
183: TSPreStep(ts);
184: TSPreStage(ts,th->stage_time);
186: if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
187: VecZeroEntries(th->Xdot);
188: if (!th->affine) {VecDuplicate(ts->vec_sol,&th->affine);}
189: TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);
190: VecScale(th->affine,(th->Theta-1.)/th->Theta);
191: }
192: if (th->extrapolate) {
193: VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);
194: } else {
195: VecCopy(ts->vec_sol,th->X);
196: }
197: SNESSolve(ts->snes,th->affine,th->X);
198: SNESGetIterationNumber(ts->snes,&its);
199: SNESGetLinearSolveIterations(ts->snes,&lits);
200: ts->snes_its += its; ts->ksp_its += lits;
201: TSPostStage(ts,th->stage_time,0,&(th->X));
202: TSGetAdapt(ts,&adapt);
203: TSAdaptCheckStage(adapt,ts,&stageok);
204: if (!stageok) {accept = PETSC_FALSE; goto reject_step;}
206: TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);
207: th->status = TS_STEP_PENDING;
208: /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
209: TSGetAdapt(ts,&adapt);
210: TSAdaptCandidatesClear(adapt);
211: TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);
212: TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);
213: if (!accept) { /* Roll back the current step */
214: ts->ptime += next_time_step; /* This will be undone in rollback */
215: th->status = TS_STEP_INCOMPLETE;
216: TSRollBack(ts);
217: goto reject_step;
218: }
220: if (ts->vec_costintegral) {
221: /* Evolve ts->vec_costintegral to compute integrals */
222: if (th->endpoint) {
223: TSAdjointComputeCostIntegrand(ts,ts->ptime,th->X0,ts->vec_costintegrand);
224: VecAXPY(ts->vec_costintegral,ts->time_step*(1.-th->Theta),ts->vec_costintegrand);
225: }
226: TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);
227: if (th->endpoint) {
228: VecAXPY(ts->vec_costintegral,ts->time_step*th->Theta,ts->vec_costintegrand);
229: }else {
230: VecAXPY(ts->vec_costintegral,ts->time_step,ts->vec_costintegrand);
231: }
232: }
234: /* ignore next_scheme for now */
235: ts->ptime += ts->time_step;
236: ts->time_step = next_time_step;
237: ts->steps++;
238: th->status = TS_STEP_COMPLETE;
239: break;
241: reject_step:
242: if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) {
243: ts->reason = TS_DIVERGED_STEP_REJECTED;
244: PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);
245: }
246: continue;
247: }
248: return(0);
249: }
253: static PetscErrorCode TSAdjointStep_Theta(TS ts)
254: {
255: TS_Theta *th = (TS_Theta*)ts->data;
256: Vec *VecDeltaLam = th->VecDeltaLam,*VecDeltaMu = th->VecDeltaMu,*VecSensiTemp = th->VecSensiTemp;
257: PetscInt nadj;
258: PetscErrorCode ierr;
259: Mat J,Jp;
260: KSP ksp;
261: PetscReal shift;
264:
265: th->status = TS_STEP_INCOMPLETE;
266: SNESGetKSP(ts->snes,&ksp);
267: TSGetIJacobian(ts,&J,&Jp,NULL,NULL);
268: th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/
270: TSPreStep(ts);
272: /* Build RHS */
273: if (ts->vec_costintegral) { /* Cost function has an integral term */
274: if (th->endpoint) {
275: TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);
276: }else {
277: TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);
278: }
279: }
280: for (nadj=0; nadj<ts->numcost; nadj++) {
281: VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);
282: VecScale(VecSensiTemp[nadj],-1./(th->Theta*ts->time_step));
283: if (ts->vec_costintegral) {
284: VecAXPY(VecSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);
285: }
286: }
287:
288: /* Build LHS */
289: shift = -1./(th->Theta*ts->time_step);
290: if (th->endpoint) {
291: TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);
292: }else {
293: TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);
294: }
295: KSPSetOperators(ksp,J,Jp);
297: /* Solve LHS X = RHS */
298: for (nadj=0; nadj<ts->numcost; nadj++) {
299: KSPSolveTranspose(ksp,VecSensiTemp[nadj],VecDeltaLam[nadj]);
300: }
302: /* Update sensitivities, and evaluate integrals if there is any */
303: if(th->endpoint && th->Theta!=1.) { /* two-stage case */
304: shift = -1./((th->Theta-1.)*ts->time_step);
305: TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);
306: if (ts->vec_costintegral) {
307: TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);
308: if (!ts->costintegralfwd) {
309: /* Evolve ts->vec_costintegral to compute integrals */
310: TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);
311: VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);
312: TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);
313: VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1.),ts->vec_costintegrand);
314: }
315: }
316: for (nadj=0; nadj<ts->numcost; nadj++) {
317: MatMultTranspose(J,VecDeltaLam[nadj],ts->vecs_sensi[nadj]);
318: if (ts->vec_costintegral) {
319: VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);
320: }
321: VecScale(ts->vecs_sensi[nadj],1./shift);
322: }
323:
324: if (ts->vecs_sensip) { /* sensitivities wrt parameters */
325: TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);
326: for (nadj=0; nadj<ts->numcost; nadj++) {
327: MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);
328: VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecDeltaMu[nadj]);
329: }
330: TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);
331: for (nadj=0; nadj<ts->numcost; nadj++) {
332: MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);
333: VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecDeltaMu[nadj]);
334: }
335: if (ts->vec_costintegral) {
336: TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);
337: for (nadj=0; nadj<ts->numcost; nadj++) {
338: VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);
339: }
340: TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);
341: for (nadj=0; nadj<ts->numcost; nadj++) {
342: VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);
343: }
344: }
345: }
346: }else { /* one-stage case */
347: shift = 0.0;
348: TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE); /* get -f_y */
349: if (ts->vec_costintegral) {
350: TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);
351: if (!ts->costintegralfwd) {
352: /* Evolve ts->vec_costintegral to compute integrals */
353: TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);
354: VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);
355: }
356: }
357: /* When th->endpoint is true and th->Theta==1 (beuler method), the Jacobian is supposed to be evaluated at ts->ptime like this:
358: if(th->endpoint) {
359: TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);
360: }
361: but ts->ptime and ts->vec_sol have the same values as th->stage_time and th->X in this case. So the code is simplified here.
362: */
363: for (nadj=0; nadj<ts->numcost; nadj++) {
364: MatMultTranspose(J,VecDeltaLam[nadj],VecSensiTemp[nadj]);
365: VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecSensiTemp[nadj]);
366: if (ts->vec_costintegral) {
367: VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);
368: }
369: }
370: if (ts->vecs_sensip) {
371: TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);
372: for (nadj=0; nadj<ts->numcost; nadj++) {
373: MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);
374: VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecDeltaMu[nadj]);
375: }
376: if (ts->vec_costintegral) {
377: TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);
378: for (nadj=0; nadj<ts->numcost; nadj++) {
379: VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);
380: }
381: }
382: }
383: }
384:
385: ts->ptime += ts->time_step;
386: ts->steps++;
387: th->status = TS_STEP_COMPLETE;
388: return(0);
389: }
393: static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
394: {
395: TS_Theta *th = (TS_Theta*)ts->data;
396: PetscReal alpha = t - ts->ptime;
400: VecCopy(ts->vec_sol,th->X);
401: if (th->endpoint) alpha *= th->Theta;
402: VecWAXPY(X,alpha,th->Xdot,th->X);
403: return(0);
404: }
406: /*------------------------------------------------------------*/
409: static PetscErrorCode TSReset_Theta(TS ts)
410: {
411: TS_Theta *th = (TS_Theta*)ts->data;
415: VecDestroy(&th->X);
416: VecDestroy(&th->Xdot);
417: VecDestroy(&th->X0);
418: VecDestroy(&th->affine);
419: VecDestroyVecs(ts->numcost,&th->VecDeltaLam);
420: VecDestroyVecs(ts->numcost,&th->VecDeltaMu);
421: VecDestroyVecs(ts->numcost,&th->VecSensiTemp);
422: return(0);
423: }
427: static PetscErrorCode TSDestroy_Theta(TS ts)
428: {
432: TSReset_Theta(ts);
433: PetscFree(ts->data);
434: PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);
435: PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);
436: PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);
437: PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);
438: return(0);
439: }
441: /*
442: This defines the nonlinear equation that is to be solved with SNES
443: G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
444: */
447: static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
448: {
449: TS_Theta *th = (TS_Theta*)ts->data;
451: Vec X0,Xdot;
452: DM dm,dmsave;
453: PetscReal shift = 1./(th->Theta*ts->time_step);
456: SNESGetDM(snes,&dm);
457: /* When using the endpoint variant, this is actually 1/Theta * Xdot */
458: TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);
459: VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);
461: /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
462: dmsave = ts->dm;
463: ts->dm = dm;
464: TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);
465: ts->dm = dmsave;
466: TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);
467: return(0);
468: }
472: static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
473: {
474: TS_Theta *th = (TS_Theta*)ts->data;
476: Vec Xdot;
477: DM dm,dmsave;
478: PetscReal shift = 1./(th->Theta*ts->time_step);
481: SNESGetDM(snes,&dm);
483: /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
484: TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);
486: dmsave = ts->dm;
487: ts->dm = dm;
488: TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);
489: ts->dm = dmsave;
490: TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);
491: return(0);
492: }
496: static PetscErrorCode TSSetUp_Theta(TS ts)
497: {
498: TS_Theta *th = (TS_Theta*)ts->data;
500: SNES snes;
501: TSAdapt adapt;
502: DM dm;
505: if (!th->X) {
506: VecDuplicate(ts->vec_sol,&th->X);
507: }
508: if (!th->Xdot) {
509: VecDuplicate(ts->vec_sol,&th->Xdot);
510: }
511: if (!th->X0) {
512: VecDuplicate(ts->vec_sol,&th->X0);
513: }
514: TSGetSNES(ts,&snes);
515: TSGetDM(ts,&dm);
516: if (dm) {
517: DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);
518: DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);
519: }
520: if (th->Theta == 0.5 && th->endpoint) th->order = 2;
521: else th->order = 1;
523: TSGetAdapt(ts,&adapt);
524: if (!th->adapt) {
525: TSAdaptSetType(adapt,TSADAPTNONE);
526: }
527: return(0);
528: }
529: /*------------------------------------------------------------*/
533: static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
534: {
535: TS_Theta *th = (TS_Theta*)ts->data;
539: VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecDeltaLam);
540: if(ts->vecs_sensip) {
541: VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecDeltaMu);
542: }
543: VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecSensiTemp);
544: return(0);
545: }
546: /*------------------------------------------------------------*/
550: static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts)
551: {
552: TS_Theta *th = (TS_Theta*)ts->data;
556: PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");
557: {
558: PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);
559: PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);
560: PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);
561: PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);
562: SNESSetFromOptions(ts->snes);
563: }
564: PetscOptionsTail();
565: return(0);
566: }
570: static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
571: {
572: TS_Theta *th = (TS_Theta*)ts->data;
573: PetscBool iascii;
577: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
578: if (iascii) {
579: PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);
580: PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");
581: }
582: if (ts->snes) {SNESView(ts->snes,viewer);}
583: return(0);
584: }
588: PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
589: {
590: TS_Theta *th = (TS_Theta*)ts->data;
593: *theta = th->Theta;
594: return(0);
595: }
599: PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
600: {
601: TS_Theta *th = (TS_Theta*)ts->data;
604: if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
605: th->Theta = theta;
606: return(0);
607: }
611: PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
612: {
613: TS_Theta *th = (TS_Theta*)ts->data;
616: *endpoint = th->endpoint;
617: return(0);
618: }
622: PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
623: {
624: TS_Theta *th = (TS_Theta*)ts->data;
627: th->endpoint = flg;
628: return(0);
629: }
631: #if defined(PETSC_HAVE_COMPLEX)
634: static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
635: {
636: PetscComplex z = xr + xi*PETSC_i,f;
637: TS_Theta *th = (TS_Theta*)ts->data;
638: const PetscReal one = 1.0;
641: f = (one + (one - th->Theta)*z)/(one - th->Theta*z);
642: *yr = PetscRealPartComplex(f);
643: *yi = PetscImaginaryPartComplex(f);
644: return(0);
645: }
646: #endif
650: static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
651: {
652: TS_Theta *th = (TS_Theta*)ts->data;
655: *ns = 1;
656: if(Y) {
657: *Y = &(th->X);
658: }
659: return(0);
660: }
662: /* ------------------------------------------------------------ */
663: /*MC
664: TSTHETA - DAE solver using the implicit Theta method
666: Level: beginner
668: Options Database:
669: -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
670: -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
671: -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
673: Notes:
674: $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
675: $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
676: $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
680: This method can be applied to DAE.
682: This method is cast as a 1-stage implicit Runge-Kutta method.
684: .vb
685: Theta | Theta
686: -------------
687: | 1
688: .ve
690: For the default Theta=0.5, this is also known as the implicit midpoint rule.
692: When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
694: .vb
695: 0 | 0 0
696: 1 | 1-Theta Theta
697: -------------------
698: | 1-Theta Theta
699: .ve
701: For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
703: To apply a diagonally implicit RK method to DAE, the stage formula
705: $ Y_i = X + h sum_j a_ij Y'_j
707: is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
709: .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
711: M*/
714: PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
715: {
716: TS_Theta *th;
720: ts->ops->reset = TSReset_Theta;
721: ts->ops->destroy = TSDestroy_Theta;
722: ts->ops->view = TSView_Theta;
723: ts->ops->setup = TSSetUp_Theta;
724: ts->ops->adjointsetup = TSAdjointSetUp_Theta;
725: ts->ops->step = TSStep_Theta;
726: ts->ops->interpolate = TSInterpolate_Theta;
727: ts->ops->evaluatestep = TSEvaluateStep_Theta;
728: ts->ops->rollback = TSRollBack_Theta;
729: ts->ops->setfromoptions = TSSetFromOptions_Theta;
730: ts->ops->snesfunction = SNESTSFormFunction_Theta;
731: ts->ops->snesjacobian = SNESTSFormJacobian_Theta;
732: #if defined(PETSC_HAVE_COMPLEX)
733: ts->ops->linearstability = TSComputeLinearStability_Theta;
734: #endif
735: ts->ops->getstages = TSGetStages_Theta;
736: ts->ops->adjointstep = TSAdjointStep_Theta;
738: PetscNewLog(ts,&th);
739: ts->data = (void*)th;
741: th->extrapolate = PETSC_FALSE;
742: th->Theta = 0.5;
743: th->ccfl = 1.0;
744: th->adapt = PETSC_FALSE;
745: PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);
746: PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);
747: PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);
748: PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);
749: return(0);
750: }
754: /*@
755: TSThetaGetTheta - Get the abscissa of the stage in (0,1].
757: Not Collective
759: Input Parameter:
760: . ts - timestepping context
762: Output Parameter:
763: . theta - stage abscissa
765: Note:
766: Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
768: Level: Advanced
770: .seealso: TSThetaSetTheta()
771: @*/
772: PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta)
773: {
779: PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));
780: return(0);
781: }
785: /*@
786: TSThetaSetTheta - Set the abscissa of the stage in (0,1].
788: Not Collective
790: Input Parameter:
791: + ts - timestepping context
792: - theta - stage abscissa
794: Options Database:
795: . -ts_theta_theta <theta>
797: Level: Intermediate
799: .seealso: TSThetaGetTheta()
800: @*/
801: PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta)
802: {
807: PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));
808: return(0);
809: }
813: /*@
814: TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
816: Not Collective
818: Input Parameter:
819: . ts - timestepping context
821: Output Parameter:
822: . endpoint - PETSC_TRUE when using the endpoint variant
824: Level: Advanced
826: .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
827: @*/
828: PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
829: {
835: PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));
836: return(0);
837: }
841: /*@
842: TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
844: Not Collective
846: Input Parameter:
847: + ts - timestepping context
848: - flg - PETSC_TRUE to use the endpoint variant
850: Options Database:
851: . -ts_theta_endpoint <flg>
853: Level: Intermediate
855: .seealso: TSTHETA, TSCN
856: @*/
857: PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
858: {
863: PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));
864: return(0);
865: }
867: /*
868: * TSBEULER and TSCN are straightforward specializations of TSTHETA.
869: * The creation functions for these specializations are below.
870: */
874: static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
875: {
879: SNESView(ts->snes,viewer);
880: return(0);
881: }
883: /*MC
884: TSBEULER - ODE solver using the implicit backward Euler method
886: Level: beginner
888: Notes:
889: TSBEULER is equivalent to TSTHETA with Theta=1.0
891: $ -ts_type theta -ts_theta_theta 1.
893: .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
895: M*/
898: PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
899: {
903: TSCreate_Theta(ts);
904: TSThetaSetTheta(ts,1.0);
905: ts->ops->view = TSView_BEuler;
906: return(0);
907: }
911: static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
912: {
916: SNESView(ts->snes,viewer);
917: return(0);
918: }
920: /*MC
921: TSCN - ODE solver using the implicit Crank-Nicolson method.
923: Level: beginner
925: Notes:
926: TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
928: $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
930: .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
932: M*/
935: PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
936: {
940: TSCreate_Theta(ts);
941: TSThetaSetTheta(ts,0.5);
942: TSThetaSetEndpoint(ts,PETSC_TRUE);
943: ts->ops->view = TSView_CN;
944: return(0);
945: }