Actual source code: theta.c

petsc-3.4.5 2014-06-29
  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 <petscsnesfas.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 TSStep_Theta(TS ts)
153: {
154:   TS_Theta            *th = (TS_Theta*)ts->data;
155:   PetscInt            its,lits,reject,next_scheme;
156:   PetscReal           next_time_step;
157:   SNESConvergedReason snesreason;
158:   PetscErrorCode      ierr;
159:   TSAdapt             adapt;
160:   PetscBool           accept;

163:   th->status = TS_STEP_INCOMPLETE;
164:   VecCopy(ts->vec_sol,th->X0);
165:   for (reject=0; reject<ts->max_reject && !ts->reason && th->status != TS_STEP_COMPLETE; reject++,ts->reject++) {
166:     PetscReal shift = 1./(th->Theta*ts->time_step);
167:     next_time_step = ts->time_step;
168:     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
169:     TSPreStep(ts);
170:     TSPreStage(ts,th->stage_time);

172:     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
173:       VecZeroEntries(th->Xdot);
174:       if (!th->affine) {VecDuplicate(ts->vec_sol,&th->affine);}
175:       TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);
176:       VecScale(th->affine,(th->Theta-1.)/th->Theta);
177:     }
178:     if (th->extrapolate) {
179:       VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);
180:     } else {
181:       VecCopy(ts->vec_sol,th->X);
182:     }
183:     SNESSolve(ts->snes,th->affine,th->X);
184:     SNESGetIterationNumber(ts->snes,&its);
185:     SNESGetLinearSolveIterations(ts->snes,&lits);
186:     SNESGetConvergedReason(ts->snes,&snesreason);
187:     ts->snes_its += its; ts->ksp_its += lits;
188:     TSGetAdapt(ts,&adapt);
189:     TSAdaptCheckStage(adapt,ts,&accept);
190:     if (!accept) continue;
191:     TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);
192:     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
193:     TSGetAdapt(ts,&adapt);
194:     TSAdaptCandidatesClear(adapt);
195:     TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);
196:     TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);

198:     if (accept) {
199:       /* ignore next_scheme for now */
200:       ts->ptime    += ts->time_step;
201:       ts->time_step = next_time_step;
202:       ts->steps++;
203:       th->status = TS_STEP_COMPLETE;
204:     } else {                    /* Roll back the current step */
205:       VecCopy(th->X0,ts->vec_sol);
206:       ts->time_step = next_time_step;
207:       th->status    = TS_STEP_INCOMPLETE;
208:     }
209:   }
210:   return(0);
211: }

215: static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
216: {
217:   TS_Theta       *th   = (TS_Theta*)ts->data;
218:   PetscReal      alpha = t - ts->ptime;

222:   VecCopy(ts->vec_sol,th->X);
223:   if (th->endpoint) alpha *= th->Theta;
224:   VecWAXPY(X,alpha,th->Xdot,th->X);
225:   return(0);
226: }

228: /*------------------------------------------------------------*/
231: static PetscErrorCode TSReset_Theta(TS ts)
232: {
233:   TS_Theta       *th = (TS_Theta*)ts->data;

237:   VecDestroy(&th->X);
238:   VecDestroy(&th->Xdot);
239:   VecDestroy(&th->X0);
240:   VecDestroy(&th->affine);
241:   return(0);
242: }

246: static PetscErrorCode TSDestroy_Theta(TS ts)
247: {

251:   TSReset_Theta(ts);
252:   PetscFree(ts->data);
253:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);
254:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);
255:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);
256:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);
257:   return(0);
258: }

260: /*
261:   This defines the nonlinear equation that is to be solved with SNES
262:   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
263: */
266: static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
267: {
268:   TS_Theta       *th = (TS_Theta*)ts->data;
270:   Vec            X0,Xdot;
271:   DM             dm,dmsave;
272:   PetscReal      shift = 1./(th->Theta*ts->time_step);

275:   SNESGetDM(snes,&dm);
276:   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
277:   TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);
278:   VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);

280:   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
281:   dmsave = ts->dm;
282:   ts->dm = dm;
283:   TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);
284:   ts->dm = dmsave;
285:   TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);
286:   return(0);
287: }

291: static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
292: {
293:   TS_Theta       *th = (TS_Theta*)ts->data;
295:   Vec            Xdot;
296:   DM             dm,dmsave;
297:   PetscReal      shift = 1./(th->Theta*ts->time_step);

300:   SNESGetDM(snes,&dm);

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

305:   dmsave = ts->dm;
306:   ts->dm = dm;
307:   TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,str,PETSC_FALSE);
308:   ts->dm = dmsave;
309:   TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);
310:   return(0);
311: }

315: static PetscErrorCode TSSetUp_Theta(TS ts)
316: {
317:   TS_Theta       *th = (TS_Theta*)ts->data;
319:   SNES           snes;
320:   DM             dm;

323:   VecDuplicate(ts->vec_sol,&th->X);
324:   VecDuplicate(ts->vec_sol,&th->Xdot);
325:   VecDuplicate(ts->vec_sol,&th->X0);
326:   TSGetSNES(ts,&snes);
327:   TSGetDM(ts,&dm);
328:   if (dm) {
329:     DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);
330:     DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);
331:   }
332:   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
333:   else th->order = 1;

335:   if (!th->adapt) {
336:     TSAdapt adapt;
337:     TSAdaptDestroy(&ts->adapt);
338:     TSGetAdapt(ts,&adapt);
339:     TSAdaptSetType(adapt,TSADAPTNONE);
340:   }
341:   return(0);
342: }
343: /*------------------------------------------------------------*/

347: static PetscErrorCode TSSetFromOptions_Theta(TS ts)
348: {
349:   TS_Theta       *th = (TS_Theta*)ts->data;

353:   PetscOptionsHead("Theta ODE solver options");
354:   {
355:     PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);
356:     PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);
357:     PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);
358:     PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);
359:     SNESSetFromOptions(ts->snes);
360:   }
361:   PetscOptionsTail();
362:   return(0);
363: }

367: static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
368: {
369:   TS_Theta       *th = (TS_Theta*)ts->data;
370:   PetscBool      iascii;

374:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
375:   if (iascii) {
376:     PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);
377:     PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");
378:   }
379:   SNESView(ts->snes,viewer);
380:   return(0);
381: }

385: PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
386: {
387:   TS_Theta *th = (TS_Theta*)ts->data;

390:   *theta = th->Theta;
391:   return(0);
392: }

396: PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
397: {
398:   TS_Theta *th = (TS_Theta*)ts->data;

401:   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
402:   th->Theta = theta;
403:   return(0);
404: }

408: PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
409: {
410:   TS_Theta *th = (TS_Theta*)ts->data;

413:   *endpoint = th->endpoint;
414:   return(0);
415: }

419: PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
420: {
421:   TS_Theta *th = (TS_Theta*)ts->data;

424:   th->endpoint = flg;
425:   return(0);
426: }

428: #if defined(PETSC_HAVE_COMPLEX)
431: static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
432: {
433:   PetscComplex z   = xr + xi*PETSC_i,f;
434:   TS_Theta     *th = (TS_Theta*)ts->data;
435:   const PetscReal one = 1.0;

438:   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
439:   *yr = PetscRealPartComplex(f);
440:   *yi = PetscImaginaryPartComplex(f);
441:   return(0);
442: }
443: #endif


446: /* ------------------------------------------------------------ */
447: /*MC
448:       TSTHETA - DAE solver using the implicit Theta method

450:    Level: beginner

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

457:    Notes:
458: $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
459: $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
460: $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)



464:    This method can be applied to DAE.

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

468: .vb
469:   Theta | Theta
470:   -------------
471:         |  1
472: .ve

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

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

478: .vb
479:   0 | 0         0
480:   1 | 1-Theta   Theta
481:   -------------------
482:     | 1-Theta   Theta
483: .ve

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

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

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

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

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

495: M*/
498: PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
499: {
500:   TS_Theta       *th;

504:   ts->ops->reset          = TSReset_Theta;
505:   ts->ops->destroy        = TSDestroy_Theta;
506:   ts->ops->view           = TSView_Theta;
507:   ts->ops->setup          = TSSetUp_Theta;
508:   ts->ops->step           = TSStep_Theta;
509:   ts->ops->interpolate    = TSInterpolate_Theta;
510:   ts->ops->evaluatestep   = TSEvaluateStep_Theta;
511:   ts->ops->setfromoptions = TSSetFromOptions_Theta;
512:   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
513:   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
514: #if defined(PETSC_HAVE_COMPLEX)
515:   ts->ops->linearstability = TSComputeLinearStability_Theta;
516: #endif

518:   PetscNewLog(ts,TS_Theta,&th);
519:   ts->data = (void*)th;

521:   th->extrapolate = PETSC_FALSE;
522:   th->Theta       = 0.5;
523:   th->ccfl        = 1.0;
524:   th->adapt       = PETSC_FALSE;
525:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);
526:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);
527:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);
528:   PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);
529:   return(0);
530: }

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

537:   Not Collective

539:   Input Parameter:
540: .  ts - timestepping context

542:   Output Parameter:
543: .  theta - stage abscissa

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

548:   Level: Advanced

550: .seealso: TSThetaSetTheta()
551: @*/
552: PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
553: {

559:   PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));
560:   return(0);
561: }

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

568:   Not Collective

570:   Input Parameter:
571: +  ts - timestepping context
572: -  theta - stage abscissa

574:   Options Database:
575: .  -ts_theta_theta <theta>

577:   Level: Intermediate

579: .seealso: TSThetaGetTheta()
580: @*/
581: PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
582: {

587:   PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));
588:   return(0);
589: }

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

596:   Not Collective

598:   Input Parameter:
599: .  ts - timestepping context

601:   Output Parameter:
602: .  endpoint - PETSC_TRUE when using the endpoint variant

604:   Level: Advanced

606: .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
607: @*/
608: PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
609: {

615:   PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));
616:   return(0);
617: }

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

624:   Not Collective

626:   Input Parameter:
627: +  ts - timestepping context
628: -  flg - PETSC_TRUE to use the endpoint variant

630:   Options Database:
631: .  -ts_theta_endpoint <flg>

633:   Level: Intermediate

635: .seealso: TSTHETA, TSCN
636: @*/
637: PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
638: {

643:   PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));
644:   return(0);
645: }

647: /*
648:  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
649:  * The creation functions for these specializations are below.
650:  */

654: static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
655: {

659:   SNESView(ts->snes,viewer);
660:   return(0);
661: }

663: /*MC
664:       TSBEULER - ODE solver using the implicit backward Euler method

666:   Level: beginner

668:   Notes:
669:   TSBEULER is equivalent to TSTHETA with Theta=1.0

671: $  -ts_type theta -ts_theta_theta 1.

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

675: M*/
678: PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
679: {

683:   TSCreate_Theta(ts);
684:   TSThetaSetTheta(ts,1.0);
685:   ts->ops->view = TSView_BEuler;
686:   return(0);
687: }

691: static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
692: {

696:   SNESView(ts->snes,viewer);
697:   return(0);
698: }

700: /*MC
701:       TSCN - ODE solver using the implicit Crank-Nicolson method.

703:   Level: beginner

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

708: $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint

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

712: M*/
715: PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
716: {

720:   TSCreate_Theta(ts);
721:   TSThetaSetTheta(ts,0.5);
722:   TSThetaSetEndpoint(ts,PETSC_TRUE);
723:   ts->ops->view = TSView_CN;
724:   return(0);
725: }