Actual source code: posindep.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: /*
  2:        Code for Timestepping with implicit backwards Euler.
  3: */
  4: #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/

  6: typedef struct {
  7:   Vec update;       /* work vector where new solution is formed */
  8:   Vec func;         /* work vector where F(t[i],u[i]) is stored */
  9:   Vec xdot;         /* work vector for time derivative of state */

 11:   /* information used for Pseudo-timestepping */

 13:   PetscErrorCode (*dt)(TS,PetscReal*,void*);              /* compute next timestep, and related context */
 14:   void *dtctx;
 15:   PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,PetscBool*);  /* verify previous timestep and related context */
 16:   void *verifyctx;

 18:   PetscReal fnorm_initial,fnorm;                   /* original and current norm of F(u) */
 19:   PetscReal fnorm_previous;

 21:   PetscReal dt_initial;                     /* initial time-step */
 22:   PetscReal dt_increment;                   /* scaling that dt is incremented each time-step */
 23:   PetscReal dt_max;                         /* maximum time step */
 24:   PetscBool increment_dt_from_initial_dt;
 25: } TS_Pseudo;

 27: /* ------------------------------------------------------------------------------*/

 31: /*@C
 32:     TSPseudoComputeTimeStep - Computes the next timestep for a currently running
 33:     pseudo-timestepping process.

 35:     Collective on TS

 37:     Input Parameter:
 38: .   ts - timestep context

 40:     Output Parameter:
 41: .   dt - newly computed timestep

 43:     Level: developer

 45:     Notes:
 46:     The routine to be called here to compute the timestep should be
 47:     set by calling TSPseudoSetTimeStep().

 49: .keywords: timestep, pseudo, compute

 51: .seealso: TSPseudoTimeStepDefault(), TSPseudoSetTimeStep()
 52: @*/
 53: PetscErrorCode  TSPseudoComputeTimeStep(TS ts,PetscReal *dt)
 54: {
 55:   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;

 59:   PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);
 60:   (*pseudo->dt)(ts,dt,pseudo->dtctx);
 61:   PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);
 62:   return(0);
 63: }


 66: /* ------------------------------------------------------------------------------*/
 69: /*@C
 70:    TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep.

 72:    Collective on TS

 74:    Input Parameters:
 75: +  ts - the timestep context
 76: .  dtctx - unused timestep context
 77: -  update - latest solution vector

 79:    Output Parameters:
 80: +  newdt - the timestep to use for the next step
 81: -  flag - flag indicating whether the last time step was acceptable

 83:    Level: advanced

 85:    Note:
 86:    This routine always returns a flag of 1, indicating an acceptable
 87:    timestep.

 89: .keywords: timestep, pseudo, default, verify

 91: .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep()
 92: @*/
 93: PetscErrorCode  TSPseudoVerifyTimeStepDefault(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool  *flag)
 94: {
 96:   *flag = PETSC_TRUE;
 97:   return(0);
 98: }


103: /*@
104:     TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable.

106:     Collective on TS

108:     Input Parameters:
109: +   ts - timestep context
110: -   update - latest solution vector

112:     Output Parameters:
113: +   dt - newly computed timestep (if it had to shrink)
114: -   flag - indicates if current timestep was ok

116:     Level: advanced

118:     Notes:
119:     The routine to be called here to compute the timestep should be
120:     set by calling TSPseudoSetVerifyTimeStep().

122: .keywords: timestep, pseudo, verify

124: .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStepDefault()
125: @*/
126: PetscErrorCode  TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool  *flag)
127: {
128:   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;

132:   if (!pseudo->verify) {*flag = PETSC_TRUE; return(0);}

134:   (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);
135:   return(0);
136: }

138: /* --------------------------------------------------------------------------------*/

142: static PetscErrorCode TSStep_Pseudo(TS ts)
143: {
144:   TS_Pseudo           *pseudo = (TS_Pseudo*)ts->data;
145:   PetscInt            its,lits,reject;
146:   PetscBool           stepok;
147:   PetscReal           next_time_step;
148:   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
149:   PetscErrorCode      ierr;

152:   if (ts->steps == 0) pseudo->dt_initial = ts->time_step;
153:   VecCopy(ts->vec_sol,pseudo->update);
154:   next_time_step = ts->time_step;
155:   TSPseudoComputeTimeStep(ts,&next_time_step);
156:   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
157:     ts->time_step = next_time_step;
158:     TSPreStep(ts);
159:     TSPreStage(ts,ts->ptime+ts->time_step);
160:     SNESSolve(ts->snes,NULL,pseudo->update);
161:     SNESGetConvergedReason(ts->snes,&snesreason);
162:     SNESGetLinearSolveIterations(ts->snes,&lits);
163:     SNESGetIterationNumber(ts->snes,&its);
164:     TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update));
165:     ts->snes_its += its; ts->ksp_its += lits;
166:     PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);
167:     pseudo->fnorm = -1;         /* The current norm is no longer valid, monitor must recompute it. */
168:     TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);
169:     if (stepok) break;
170:   }
171:   if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
172:     ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
173:     PetscInfo2(ts,"step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
174:     return(0);
175:   }
176:   if (reject >= ts->max_reject) {
177:     ts->reason = TS_DIVERGED_STEP_REJECTED;
178:     PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);
179:     return(0);
180:   }
181:   VecCopy(pseudo->update,ts->vec_sol);
182:   ts->ptime += ts->time_step;
183:   ts->time_step = next_time_step;
184:   ts->steps++;
185:   return(0);
186: }

188: /*------------------------------------------------------------*/
191: static PetscErrorCode TSReset_Pseudo(TS ts)
192: {
193:   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;

197:   VecDestroy(&pseudo->update);
198:   VecDestroy(&pseudo->func);
199:   VecDestroy(&pseudo->xdot);
200:   return(0);
201: }

205: static PetscErrorCode TSDestroy_Pseudo(TS ts)
206: {

210:   TSReset_Pseudo(ts);
211:   PetscFree(ts->data);
212:   PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);
213:   PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);
214:   PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);
215:   PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);
216:   PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);
217:   return(0);
218: }

220: /*------------------------------------------------------------*/

224: /*
225:     Compute Xdot = (X^{n+1}-X^n)/dt) = 0
226: */
227: static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot)
228: {
229:   TS_Pseudo         *pseudo = (TS_Pseudo*)ts->data;
230:   const PetscScalar mdt     = 1.0/ts->time_step,*xnp1,*xn;
231:   PetscScalar       *xdot;
232:   PetscErrorCode    ierr;
233:   PetscInt          i,n;

236:   VecGetArrayRead(ts->vec_sol,&xn);
237:   VecGetArrayRead(X,&xnp1);
238:   VecGetArray(pseudo->xdot,&xdot);
239:   VecGetLocalSize(X,&n);
240:   for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]);
241:   VecRestoreArrayRead(ts->vec_sol,&xn);
242:   VecRestoreArrayRead(X,&xnp1);
243:   VecRestoreArray(pseudo->xdot,&xdot);
244:   *Xdot = pseudo->xdot;
245:   return(0);
246: }

250: /*
251:     The transient residual is

253:         F(U^{n+1},(U^{n+1}-U^n)/dt) = 0

255:     or for ODE,

257:         (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0

259:     This is the function that must be evaluated for transient simulation and for
260:     finite difference Jacobians.  On the first Newton step, this algorithm uses
261:     a guess of U^{n+1} = U^n in which case the transient term vanishes and the
262:     residual is actually the steady state residual.  Pseudotransient
263:     continuation as described in the literature is a linearly implicit
264:     algorithm, it only takes this one Newton step with the steady state
265:     residual, and then advances to the next time step.
266: */
267: static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts)
268: {
269:   Vec            Xdot;

273:   TSPseudoGetXdot(ts,X,&Xdot);
274:   TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);
275:   return(0);
276: }

280: /*
281:    This constructs the Jacobian needed for SNES.  For DAE, this is

283:        dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot

285:     and for ODE:

287:        J = I/dt - J_{Frhs}   where J_{Frhs} is the given Jacobian of Frhs.
288: */
289: static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts)
290: {
291:   Vec            Xdot;

295:   TSPseudoGetXdot(ts,X,&Xdot);
296:   TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE);
297:   return(0);
298: }


303: static PetscErrorCode TSSetUp_Pseudo(TS ts)
304: {
305:   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;

309:   VecDuplicate(ts->vec_sol,&pseudo->update);
310:   VecDuplicate(ts->vec_sol,&pseudo->func);
311:   VecDuplicate(ts->vec_sol,&pseudo->xdot);
312:   return(0);
313: }
314: /*------------------------------------------------------------*/

318: PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
319: {
320:   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
322:   PetscViewer    viewer = (PetscViewer) dummy;

325:   if (!viewer) {
326:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);
327:   }
328:   if (pseudo->fnorm < 0) {      /* The last computed norm is stale, recompute */
329:     VecZeroEntries(pseudo->xdot);
330:     TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);
331:     VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);
332:   }
333:   PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);
334:   PetscViewerASCIIPrintf(viewer,"TS %D dt %g time %g fnorm %g\n",step,(double)ts->time_step,(double)ptime,(double)pseudo->fnorm);
335:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);
336:   return(0);
337: }

341: static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptions *PetscOptionsObject,TS ts)
342: {
343:   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
345:   PetscBool      flg = PETSC_FALSE;
346:   PetscViewer    viewer;

349:   PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options");
350:   PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,NULL);
351:   if (flg) {
352:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);
353:     TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
354:   }
355:   flg  = PETSC_FALSE;
356:   PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);
357:   if (flg) {
358:     TSPseudoIncrementDtFromInitialDt(ts);
359:   }
360:   PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL);
361:   PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL);

363:   SNESSetFromOptions(ts->snes);
364:   PetscOptionsTail();
365:   return(0);
366: }

370: static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer)
371: {

375:   SNESView(ts->snes,viewer);
376:   return(0);
377: }

379: /* ----------------------------------------------------------------------------- */
382: /*@C
383:    TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the
384:    last timestep.

386:    Logically Collective on TS

388:    Input Parameters:
389: +  ts - timestep context
390: .  dt - user-defined function to verify timestep
391: -  ctx - [optional] user-defined context for private data
392:          for the timestep verification routine (may be NULL)

394:    Level: advanced

396:    Calling sequence of func:
397: .  func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool  *flag);

399: .  update - latest solution vector
400: .  ctx - [optional] timestep context
401: .  newdt - the timestep to use for the next step
402: .  flag - flag indicating whether the last time step was acceptable

404:    Notes:
405:    The routine set here will be called by TSPseudoVerifyTimeStep()
406:    during the timestepping process.

408: .keywords: timestep, pseudo, set, verify

410: .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep()
411: @*/
412: PetscErrorCode  TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx)
413: {

418:   PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));
419:   return(0);
420: }

424: /*@
425:     TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to
426:     dt when using the TSPseudoTimeStepDefault() routine.

428:    Logically Collective on TS

430:     Input Parameters:
431: +   ts - the timestep context
432: -   inc - the scaling factor >= 1.0

434:     Options Database Key:
435: .    -ts_pseudo_increment <increment>

437:     Level: advanced

439: .keywords: timestep, pseudo, set, increment

441: .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
442: @*/
443: PetscErrorCode  TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc)
444: {

450:   PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));
451:   return(0);
452: }

456: /*@
457:     TSPseudoSetMaxTimeStep - Sets the maximum time step
458:     when using the TSPseudoTimeStepDefault() routine.

460:    Logically Collective on TS

462:     Input Parameters:
463: +   ts - the timestep context
464: -   maxdt - the maximum time step, use a non-positive value to deactivate

466:     Options Database Key:
467: .    -ts_pseudo_max_dt <increment>

469:     Level: advanced

471: .keywords: timestep, pseudo, set

473: .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
474: @*/
475: PetscErrorCode  TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt)
476: {

482:   PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));
483:   return(0);
484: }

488: /*@
489:     TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep
490:     is computed via the formula
491: $         dt = initial_dt*initial_fnorm/current_fnorm
492:       rather than the default update,
493: $         dt = current_dt*previous_fnorm/current_fnorm.

495:    Logically Collective on TS

497:     Input Parameter:
498: .   ts - the timestep context

500:     Options Database Key:
501: .    -ts_pseudo_increment_dt_from_initial_dt

503:     Level: advanced

505: .keywords: timestep, pseudo, set, increment

507: .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault()
508: @*/
509: PetscErrorCode  TSPseudoIncrementDtFromInitialDt(TS ts)
510: {

515:   PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));
516:   return(0);
517: }


522: /*@C
523:    TSPseudoSetTimeStep - Sets the user-defined routine to be
524:    called at each pseudo-timestep to update the timestep.

526:    Logically Collective on TS

528:    Input Parameters:
529: +  ts - timestep context
530: .  dt - function to compute timestep
531: -  ctx - [optional] user-defined context for private data
532:          required by the function (may be NULL)

534:    Level: intermediate

536:    Calling sequence of func:
537: .  func (TS ts,PetscReal *newdt,void *ctx);

539: .  newdt - the newly computed timestep
540: .  ctx - [optional] timestep context

542:    Notes:
543:    The routine set here will be called by TSPseudoComputeTimeStep()
544:    during the timestepping process.
545:    If not set then TSPseudoTimeStepDefault() is automatically used

547: .keywords: timestep, pseudo, set

549: .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep()
550: @*/
551: PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx)
552: {

557:   PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));
558:   return(0);
559: }

561: /* ----------------------------------------------------------------------------- */

563: typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*);  /* force argument to next function to not be extern C*/
566: PetscErrorCode  TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx)
567: {
568:   TS_Pseudo *pseudo;

571:   pseudo            = (TS_Pseudo*)ts->data;
572:   pseudo->verify    = dt;
573:   pseudo->verifyctx = ctx;
574:   return(0);
575: }

579: PetscErrorCode  TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc)
580: {
581:   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;

584:   pseudo->dt_increment = inc;
585:   return(0);
586: }

590: PetscErrorCode  TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
591: {
592:   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;

595:   pseudo->dt_max = maxdt;
596:   return(0);
597: }

601: PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
602: {
603:   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;

606:   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
607:   return(0);
608: }

610: typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
613: PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx)
614: {
615:   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;

618:   pseudo->dt    = dt;
619:   pseudo->dtctx = ctx;
620:   return(0);
621: }

623: /* ----------------------------------------------------------------------------- */
624: /*MC
625:       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping

627:   This method solves equations of the form

629: $    F(X,Xdot) = 0

631:   for steady state using the iteration

633: $    [G'] S = -F(X,0)
634: $    X += S

636:   where

638: $    G(Y) = F(Y,(Y-X)/dt)

640:   This is linearly-implicit Euler with the residual always evaluated "at steady
641:   state".  See note below.

643:   Options database keys:
644: +  -ts_pseudo_increment <real> - ratio of increase dt
645: -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt

647:   Level: beginner

649:   References:
650:   Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003.
651:   C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.

653:   Notes:
654:   The residual computed by this method includes the transient term (Xdot is computed instead of
655:   always being zero), but since the prediction from the last step is always the solution from the
656:   last step, on the first Newton iteration we have

658: $  Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0

660:   Therefore, the linear system solved by the first Newton iteration is equivalent to the one
661:   described above and in the papers.  If the user chooses to perform multiple Newton iterations, the
662:   algorithm is no longer the one described in the referenced papers.

664: .seealso:  TSCreate(), TS, TSSetType()

666: M*/
669: PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts)
670: {
671:   TS_Pseudo      *pseudo;
673:   SNES           snes;
674:   SNESType       stype;

677:   ts->ops->reset   = TSReset_Pseudo;
678:   ts->ops->destroy = TSDestroy_Pseudo;
679:   ts->ops->view    = TSView_Pseudo;

681:   ts->ops->setup          = TSSetUp_Pseudo;
682:   ts->ops->step           = TSStep_Pseudo;
683:   ts->ops->setfromoptions = TSSetFromOptions_Pseudo;
684:   ts->ops->snesfunction   = SNESTSFormFunction_Pseudo;
685:   ts->ops->snesjacobian   = SNESTSFormJacobian_Pseudo;

687:   TSGetSNES(ts,&snes);
688:   SNESGetType(snes,&stype);
689:   if (!stype) {SNESSetType(snes,SNESKSPONLY);}

691:   PetscNewLog(ts,&pseudo);
692:   ts->data = (void*)pseudo;

694:   pseudo->dt_increment                 = 1.1;
695:   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
696:   pseudo->dt                           = TSPseudoTimeStepDefault;
697:   pseudo->fnorm                        = -1;

699:   PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);
700:   PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);
701:   PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);
702:   PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);
703:   PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);
704:   return(0);
705: }

709: /*@C
710:    TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping.
711:    Use with TSPseudoSetTimeStep().

713:    Collective on TS

715:    Input Parameters:
716: .  ts - the timestep context
717: .  dtctx - unused timestep context

719:    Output Parameter:
720: .  newdt - the timestep to use for the next step

722:    Level: advanced

724: .keywords: timestep, pseudo, default

726: .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
727: @*/
728: PetscErrorCode  TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx)
729: {
730:   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
731:   PetscReal      inc     = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;

735:   VecZeroEntries(pseudo->xdot);
736:   TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);
737:   VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);
738:   if (pseudo->fnorm_initial == 0.0) {
739:     /* first time through so compute initial function norm */
740:     pseudo->fnorm_initial = pseudo->fnorm;
741:     fnorm_previous        = pseudo->fnorm;
742:   }
743:   if (pseudo->fnorm == 0.0)                      *newdt = 1.e12*inc*ts->time_step;
744:   else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
745:   else                                           *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
746:   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
747:   pseudo->fnorm_previous = pseudo->fnorm;
748:   return(0);
749: }