Actual source code: theta.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: /*
  2:   Code for timestepping with implicit Theta method
  3: */
  4: #define PETSC_DESIRE_COMPLEX
  5: #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
  6: #include <petscsnes.h>
  7: #include <petscdm.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:   PetscBool    extrapolate;
 14:   PetscBool    endpoint;
 15:   PetscReal    Theta;
 16:   PetscReal    stage_time;
 17:   TSStepStatus status;
 18:   char         *name;
 19:   PetscInt     order;
 20:   PetscReal    ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
 21:   PetscBool    adapt;  /* use time-step adaptivity ? */
 22: } TS_Theta;

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

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


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

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

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

 72:   return(0);
 73: }

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

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

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

101:   return(0);
102: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

321:   SNESGetDM(snes,&dm);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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


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

470:    Level: beginner

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

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



484:    This method can be applied to DAE.

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

558:   Not Collective

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

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

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

569:   Level: Advanced

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

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

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

589:   Not Collective

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

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

598:   Level: Intermediate

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

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

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

617:   Not Collective

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

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

625:   Level: Advanced

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

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

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

645:   Not Collective

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

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

654:   Level: Intermediate

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

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

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

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

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

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

687:   Level: beginner

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

692: $  -ts_type theta -ts_theta_theta 1.

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

696: M*/
699: PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
700: {

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

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

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

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

724:   Level: beginner

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

729: $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint

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

733: M*/
736: PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
737: {

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