Actual source code: theta.c

petsc-master 2014-12-27
Report Typos and Errors
  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>

  8: typedef struct {
  9:   Vec          X,Xdot;                   /* Storage for one stage */
 10:   Vec          X0;                       /* work vector to store X0 */
 11:   Vec          affine;                   /* Affine vector needed for residual at beginning of step */
 12:   PetscBool    extrapolate;
 13:   PetscBool    endpoint;
 14:   PetscReal    Theta;
 15:   PetscReal    stage_time;
 16:   TSStepStatus status;
 17:   char         *name;
 18:   PetscInt     order;
 19:   PetscReal    ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
 20:   PetscBool    adapt;  /* use time-step adaptivity ? */
 21: } TS_Theta;

 25: static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
 26: {
 27:   TS_Theta       *th = (TS_Theta*)ts->data;

 31:   if (X0) {
 32:     if (dm && dm != ts->dm) {
 33:       DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);
 34:     } else *X0 = ts->vec_sol;
 35:   }
 36:   if (Xdot) {
 37:     if (dm && dm != ts->dm) {
 38:       DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);
 39:     } else *Xdot = th->Xdot;
 40:   }
 41:   return(0);
 42: }


 47: static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
 48: {

 52:   if (X0) {
 53:     if (dm && dm != ts->dm) {
 54:       DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);
 55:     }
 56:   }
 57:   if (Xdot) {
 58:     if (dm && dm != ts->dm) {
 59:       DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);
 60:     }
 61:   }
 62:   return(0);
 63: }

 67: static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
 68: {

 71:   return(0);
 72: }

 76: static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
 77: {
 78:   TS             ts = (TS)ctx;
 80:   Vec            X0,Xdot,X0_c,Xdot_c;

 83:   TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);
 84:   TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);
 85:   MatRestrict(restrct,X0,X0_c);
 86:   MatRestrict(restrct,Xdot,Xdot_c);
 87:   VecPointwiseMult(X0_c,rscale,X0_c);
 88:   VecPointwiseMult(Xdot_c,rscale,Xdot_c);
 89:   TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);
 90:   TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);
 91:   return(0);
 92: }

 96: static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
 97: {

100:   return(0);
101: }

105: static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
106: {
107:   TS             ts = (TS)ctx;
109:   Vec            X0,Xdot,X0_sub,Xdot_sub;

112:   TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);
113:   TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);

115:   VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);
116:   VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);

118:   VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);
119:   VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);

121:   TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);
122:   TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);
123:   return(0);
124: }

128: static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done)
129: {
131:   TS_Theta       *th = (TS_Theta*)ts->data;

134:   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");
135:   if (order == th->order) {
136:     if (th->endpoint) {
137:       VecCopy(th->X,U);
138:     } else {
139:       PetscReal shift = 1./(th->Theta*ts->time_step);
140:       VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);
141:       VecAXPY(U,ts->time_step,th->Xdot);
142:     }
143:   } else if (order == th->order-1 && order) {
144:     VecWAXPY(U,ts->time_step,th->Xdot,th->X0);
145:   }
146:   return(0);
147: }

151: static PetscErrorCode TSRollBack_Theta(TS ts)
152: {
153:   TS_Theta       *th = (TS_Theta*)ts->data;

157:   VecCopy(th->X0,ts->vec_sol);
158:   th->status    = TS_STEP_INCOMPLETE;
159:   return(0);
160: }

164: static PetscErrorCode TSStep_Theta(TS ts)
165: {
166:   TS_Theta       *th = (TS_Theta*)ts->data;
167:   PetscInt       its,lits,reject,next_scheme;
168:   PetscReal      next_time_step;
169:   TSAdapt        adapt;
170:   PetscBool      stageok,accept = PETSC_TRUE;

174:   th->status = TS_STEP_INCOMPLETE;
175:   VecCopy(ts->vec_sol,th->X0);
176:   for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) {
177:     PetscReal shift = 1./(th->Theta*ts->time_step);
178:     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
179:     TSPreStep(ts);
180:     TSPreStage(ts,th->stage_time);

182:     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
183:       VecZeroEntries(th->Xdot);
184:       if (!th->affine) {VecDuplicate(ts->vec_sol,&th->affine);}
185:       TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);
186:       VecScale(th->affine,(th->Theta-1.)/th->Theta);
187:     }
188:     if (th->extrapolate) {
189:       VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);
190:     } else {
191:       VecCopy(ts->vec_sol,th->X);
192:     }
193:     SNESSolve(ts->snes,th->affine,th->X);
194:     SNESGetIterationNumber(ts->snes,&its);
195:     SNESGetLinearSolveIterations(ts->snes,&lits);
196:     ts->snes_its += its; ts->ksp_its += lits;
197:     TSPostStage(ts,th->stage_time,0,&(th->X));
198:     TSGetAdapt(ts,&adapt);
199:     TSAdaptCheckStage(adapt,ts,&stageok);
200:     if (!stageok) {accept = PETSC_FALSE; goto reject_step;}

202:     TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);
203:     th->status = TS_STEP_PENDING;
204:     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
205:     TSGetAdapt(ts,&adapt);
206:     TSAdaptCandidatesClear(adapt);
207:     TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);
208:     TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);
209:     if (!accept) {           /* Roll back the current step */
210:       ts->ptime += next_time_step; /* This will be undone in rollback */
211:       th->status = TS_STEP_INCOMPLETE;
212:       TSRollBack(ts);
213:       goto reject_step;
214:     }

216:     /* ignore next_scheme for now */
217:     ts->ptime    += ts->time_step;
218:     ts->time_step = next_time_step;
219:     ts->steps++;
220:     th->status = TS_STEP_COMPLETE;
221:     break;

223: reject_step:
224:     if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) {
225:       ts->reason = TS_DIVERGED_STEP_REJECTED;
226:       PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);
227:     }
228:     continue;
229:   }
230:   return(0);
231: }

235: static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
236: {
237:   TS_Theta       *th   = (TS_Theta*)ts->data;
238:   PetscReal      alpha = t - ts->ptime;

242:   VecCopy(ts->vec_sol,th->X);
243:   if (th->endpoint) alpha *= th->Theta;
244:   VecWAXPY(X,alpha,th->Xdot,th->X);
245:   return(0);
246: }

248: /*------------------------------------------------------------*/
251: static PetscErrorCode TSReset_Theta(TS ts)
252: {
253:   TS_Theta       *th = (TS_Theta*)ts->data;

257:   VecDestroy(&th->X);
258:   VecDestroy(&th->Xdot);
259:   VecDestroy(&th->X0);
260:   VecDestroy(&th->affine);
261:   return(0);
262: }

266: static PetscErrorCode TSDestroy_Theta(TS ts)
267: {

271:   TSReset_Theta(ts);
272:   PetscFree(ts->data);
273:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);
274:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);
275:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);
276:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);
277:   return(0);
278: }

280: /*
281:   This defines the nonlinear equation that is to be solved with SNES
282:   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
283: */
286: static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
287: {
288:   TS_Theta       *th = (TS_Theta*)ts->data;
290:   Vec            X0,Xdot;
291:   DM             dm,dmsave;
292:   PetscReal      shift = 1./(th->Theta*ts->time_step);

295:   SNESGetDM(snes,&dm);
296:   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
297:   TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);
298:   VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);

300:   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
301:   dmsave = ts->dm;
302:   ts->dm = dm;
303:   TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);
304:   ts->dm = dmsave;
305:   TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);
306:   return(0);
307: }

311: static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
312: {
313:   TS_Theta       *th = (TS_Theta*)ts->data;
315:   Vec            Xdot;
316:   DM             dm,dmsave;
317:   PetscReal      shift = 1./(th->Theta*ts->time_step);

320:   SNESGetDM(snes,&dm);

322:   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
323:   TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);

325:   dmsave = ts->dm;
326:   ts->dm = dm;
327:   TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);
328:   ts->dm = dmsave;
329:   TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);
330:   return(0);
331: }

335: static PetscErrorCode TSSetUp_Theta(TS ts)
336: {
337:   TS_Theta       *th = (TS_Theta*)ts->data;
339:   SNES           snes;
340:   TSAdapt        adapt;
341:   DM             dm;

344:   VecDuplicate(ts->vec_sol,&th->X);
345:   VecDuplicate(ts->vec_sol,&th->Xdot);
346:   VecDuplicate(ts->vec_sol,&th->X0);
347:   TSGetSNES(ts,&snes);
348:   TSGetDM(ts,&dm);
349:   if (dm) {
350:     DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);
351:     DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);
352:   }
353:   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
354:   else th->order = 1;

356:   TSGetAdapt(ts,&adapt);
357:   if (!th->adapt) {
358:     TSAdaptSetType(adapt,TSADAPTNONE);
359:   }
360:   return(0);
361: }
362: /*------------------------------------------------------------*/

366: static PetscErrorCode TSSetFromOptions_Theta(TS ts)
367: {
368:   TS_Theta       *th = (TS_Theta*)ts->data;

372:   PetscOptionsHead("Theta ODE solver options");
373:   {
374:     PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);
375:     PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);
376:     PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);
377:     PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);
378:     SNESSetFromOptions(ts->snes);
379:   }
380:   PetscOptionsTail();
381:   return(0);
382: }

386: static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
387: {
388:   TS_Theta       *th = (TS_Theta*)ts->data;
389:   PetscBool      iascii;

393:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
394:   if (iascii) {
395:     PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);
396:     PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");
397:   }
398:   if (ts->snes) {SNESView(ts->snes,viewer);}
399:   return(0);
400: }

404: PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
405: {
406:   TS_Theta *th = (TS_Theta*)ts->data;

409:   *theta = th->Theta;
410:   return(0);
411: }

415: PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
416: {
417:   TS_Theta *th = (TS_Theta*)ts->data;

420:   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
421:   th->Theta = theta;
422:   return(0);
423: }

427: PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
428: {
429:   TS_Theta *th = (TS_Theta*)ts->data;

432:   *endpoint = th->endpoint;
433:   return(0);
434: }

438: PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
439: {
440:   TS_Theta *th = (TS_Theta*)ts->data;

443:   th->endpoint = flg;
444:   return(0);
445: }

447: #if defined(PETSC_HAVE_COMPLEX)
450: static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
451: {
452:   PetscComplex z   = xr + xi*PETSC_i,f;
453:   TS_Theta     *th = (TS_Theta*)ts->data;
454:   const PetscReal one = 1.0;

457:   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
458:   *yr = PetscRealPartComplex(f);
459:   *yi = PetscImaginaryPartComplex(f);
460:   return(0);
461: }
462: #endif


465: /* ------------------------------------------------------------ */
466: /*MC
467:       TSTHETA - DAE solver using the implicit Theta method

469:    Level: beginner

471:    Options Database:
472:       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
473:       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
474:       -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method

476:    Notes:
477: $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
478: $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
479: $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)



483:    This method can be applied to DAE.

485:    This method is cast as a 1-stage implicit Runge-Kutta method.

487: .vb
488:   Theta | Theta
489:   -------------
490:         |  1
491: .ve

493:    For the default Theta=0.5, this is also known as the implicit midpoint rule.

495:    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:

497: .vb
498:   0 | 0         0
499:   1 | 1-Theta   Theta
500:   -------------------
501:     | 1-Theta   Theta
502: .ve

504:    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).

506:    To apply a diagonally implicit RK method to DAE, the stage formula

508: $  Y_i = X + h sum_j a_ij Y'_j

510:    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)

512: .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()

514: M*/
517: PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
518: {
519:   TS_Theta       *th;

523:   ts->ops->reset          = TSReset_Theta;
524:   ts->ops->destroy        = TSDestroy_Theta;
525:   ts->ops->view           = TSView_Theta;
526:   ts->ops->setup          = TSSetUp_Theta;
527:   ts->ops->step           = TSStep_Theta;
528:   ts->ops->interpolate    = TSInterpolate_Theta;
529:   ts->ops->evaluatestep   = TSEvaluateStep_Theta;
530:   ts->ops->rollback       = TSRollBack_Theta;
531:   ts->ops->setfromoptions = TSSetFromOptions_Theta;
532:   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
533:   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
534: #if defined(PETSC_HAVE_COMPLEX)
535:   ts->ops->linearstability = TSComputeLinearStability_Theta;
536: #endif

538:   PetscNewLog(ts,&th);
539:   ts->data = (void*)th;

541:   th->extrapolate = PETSC_FALSE;
542:   th->Theta       = 0.5;
543:   th->ccfl        = 1.0;
544:   th->adapt       = PETSC_FALSE;
545:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);
546:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);
547:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);
548:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);
549:   return(0);
550: }

554: /*@
555:   TSThetaGetTheta - Get the abscissa of the stage in (0,1].

557:   Not Collective

559:   Input Parameter:
560: .  ts - timestepping context

562:   Output Parameter:
563: .  theta - stage abscissa

565:   Note:
566:   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.

568:   Level: Advanced

570: .seealso: TSThetaSetTheta()
571: @*/
572: PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
573: {

579:   PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));
580:   return(0);
581: }

585: /*@
586:   TSThetaSetTheta - Set the abscissa of the stage in (0,1].

588:   Not Collective

590:   Input Parameter:
591: +  ts - timestepping context
592: -  theta - stage abscissa

594:   Options Database:
595: .  -ts_theta_theta <theta>

597:   Level: Intermediate

599: .seealso: TSThetaGetTheta()
600: @*/
601: PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
602: {

607:   PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));
608:   return(0);
609: }

613: /*@
614:   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).

616:   Not Collective

618:   Input Parameter:
619: .  ts - timestepping context

621:   Output Parameter:
622: .  endpoint - PETSC_TRUE when using the endpoint variant

624:   Level: Advanced

626: .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
627: @*/
628: PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
629: {

635:   PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));
636:   return(0);
637: }

641: /*@
642:   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).

644:   Not Collective

646:   Input Parameter:
647: +  ts - timestepping context
648: -  flg - PETSC_TRUE to use the endpoint variant

650:   Options Database:
651: .  -ts_theta_endpoint <flg>

653:   Level: Intermediate

655: .seealso: TSTHETA, TSCN
656: @*/
657: PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
658: {

663:   PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));
664:   return(0);
665: }

667: /*
668:  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
669:  * The creation functions for these specializations are below.
670:  */

674: static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
675: {

679:   SNESView(ts->snes,viewer);
680:   return(0);
681: }

683: /*MC
684:       TSBEULER - ODE solver using the implicit backward Euler method

686:   Level: beginner

688:   Notes:
689:   TSBEULER is equivalent to TSTHETA with Theta=1.0

691: $  -ts_type theta -ts_theta_theta 1.

693: .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA

695: M*/
698: PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
699: {

703:   TSCreate_Theta(ts);
704:   TSThetaSetTheta(ts,1.0);
705:   ts->ops->view = TSView_BEuler;
706:   return(0);
707: }

711: static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
712: {

716:   SNESView(ts->snes,viewer);
717:   return(0);
718: }

720: /*MC
721:       TSCN - ODE solver using the implicit Crank-Nicolson method.

723:   Level: beginner

725:   Notes:
726:   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.

728: $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint

730: .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA

732: M*/
735: PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
736: {

740:   TSCreate_Theta(ts);
741:   TSThetaSetTheta(ts,0.5);
742:   TSThetaSetEndpoint(ts,PETSC_TRUE);
743:   ts->ops->view = TSView_CN;
744:   return(0);
745: }