Actual source code: posindep.c

petsc-3.3-p1 2012-06-15
  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: /*@
 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: advanced

 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: TSPseudoDefaultTimeStep(), 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:    TSPseudoDefaultVerifyTimeStep - 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  TSPseudoDefaultVerifyTimeStep(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(), TSPseudoDefaultVerifyTimeStep()
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);

136:   return(0);
137: }

139: /* --------------------------------------------------------------------------------*/

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

153:   if (ts->steps == 0) {
154:     pseudo->dt_initial = ts->time_step;
155:   }
156:   VecCopy(ts->vec_sol,pseudo->update);
157:   next_time_step = ts->time_step;
158:   TSPseudoComputeTimeStep(ts,&next_time_step);
159:   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
160:     ts->time_step = next_time_step;
161:     SNESSolve(ts->snes,PETSC_NULL,pseudo->update);
162:     SNESGetConvergedReason(ts->snes,&snesreason);
163:     SNESGetLinearSolveIterations(ts->snes,&lits);
164:     SNESGetIterationNumber(ts->snes,&its);
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:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C","",PETSC_NULL);
213:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C","",PETSC_NULL);
214:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetMaxTimeStep_C","",PETSC_NULL);
215:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","",PETSC_NULL);
216:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","",PETSC_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;
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++) {
241:     xdot[i] = mdt*(xnp1[i] - xn[i]);
242:   }
243:   VecRestoreArrayRead(ts->vec_sol,&xn);
244:   VecRestoreArrayRead(X,&xnp1);
245:   VecRestoreArray(pseudo->xdot,&xdot);
246:   *Xdot = pseudo->xdot;
247:   return(0);
248: }

252: /*
253:     The transient residual is

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

257:     or for ODE,

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

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

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

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

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

287:     and for ODE:

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

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


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

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

320: PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
321: {
322:   TS_Pseudo        *pseudo = (TS_Pseudo*)ts->data;
323:   PetscErrorCode   ierr;
324:   PetscViewer      viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);

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

340: static PetscErrorCode TSSetFromOptions_Pseudo(TS ts)
341: {
342:   TS_Pseudo       *pseudo = (TS_Pseudo*)ts->data;
343:   PetscErrorCode  ierr;
344:   PetscBool       flg = PETSC_FALSE;
345:   PetscViewer     viewer;

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

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

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

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

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

385:    Logically Collective on TS

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

393:    Level: advanced

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

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

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

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

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

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

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

427:    Logically Collective on TS

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

433:     Options Database Key:
434: $    -ts_pseudo_increment <increment>

436:     Level: advanced

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

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

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

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

459:    Logically Collective on TS

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

465:     Options Database Key:
466: $    -ts_pseudo_max_dt <increment>

468:     Level: advanced

470: .keywords: timestep, pseudo, set

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

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

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

494:    Logically Collective on TS

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

499:     Options Database Key:
500: $    -ts_pseudo_increment_dt_from_initial_dt

502:     Level: advanced

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

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

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


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

525:    Logically Collective on TS

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

533:    Level: intermediate

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

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

541:    Notes:
542:    The routine set here will be called by TSPseudoComputeTimeStep()
543:    during the timestepping process.

545: .keywords: timestep, pseudo, set

547: .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep()
548: @*/
549: PetscErrorCode  TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx)
550: {

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

559: /* ----------------------------------------------------------------------------- */

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

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

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

585:   pseudo->dt_increment = inc;
586:   return(0);
587: }
588: EXTERN_C_END

590: EXTERN_C_BEGIN
593: PetscErrorCode  TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt)
594: {
595:   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;

598:   pseudo->dt_max = maxdt;
599:   return(0);
600: }
601: EXTERN_C_END

603: EXTERN_C_BEGIN
606: PetscErrorCode  TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts)
607: {
608:   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;

611:   pseudo->increment_dt_from_initial_dt = PETSC_TRUE;
612:   return(0);
613: }
614: EXTERN_C_END

616: typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/
617: EXTERN_C_BEGIN
620: PetscErrorCode  TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx)
621: {
622:   TS_Pseudo *pseudo = (TS_Pseudo*)ts->data;

625:   pseudo->dt      = dt;
626:   pseudo->dtctx   = ctx;
627:   return(0);
628: }
629: EXTERN_C_END

631: /* ----------------------------------------------------------------------------- */
632: /*MC
633:       TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping

635:   This method solves equations of the form

637: $    F(X,Xdot) = 0

639:   for steady state using the iteration

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

644:   where

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

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

651:   Options database keys:
652: +  -ts_pseudo_increment <real> - ratio of increase dt
653: -  -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt

655:   Level: beginner

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

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

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

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

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

674: M*/
675: EXTERN_C_BEGIN
678: PetscErrorCode  TSCreate_Pseudo(TS ts)
679: {
680:   TS_Pseudo      *pseudo;
682:   SNES           snes;
683:   const SNESType stype;

686:   ts->ops->reset           = TSReset_Pseudo;
687:   ts->ops->destroy         = TSDestroy_Pseudo;
688:   ts->ops->view            = TSView_Pseudo;

690:   ts->ops->setup           = TSSetUp_Pseudo;
691:   ts->ops->step            = TSStep_Pseudo;
692:   ts->ops->setfromoptions  = TSSetFromOptions_Pseudo;
693:   ts->ops->snesfunction    = SNESTSFormFunction_Pseudo;
694:   ts->ops->snesjacobian    = SNESTSFormJacobian_Pseudo;

696:   TSGetSNES(ts,&snes);
697:   SNESGetType(snes,&stype);
698:   if (!stype) {SNESSetType(snes,SNESKSPONLY);}

700:   PetscNewLog(ts,TS_Pseudo,&pseudo);
701:   ts->data = (void*)pseudo;

703:   pseudo->dt_increment                 = 1.1;
704:   pseudo->increment_dt_from_initial_dt = PETSC_FALSE;
705:   pseudo->dt                           = TSPseudoDefaultTimeStep;
706:   pseudo->fnorm                        = -1;

708:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",
709:                     "TSPseudoSetVerifyTimeStep_Pseudo",
710:                      TSPseudoSetVerifyTimeStep_Pseudo);
711:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",
712:                     "TSPseudoSetTimeStepIncrement_Pseudo",
713:                      TSPseudoSetTimeStepIncrement_Pseudo);
714:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",
715:                     "TSPseudoSetMaxTimeStep_Pseudo",
716:                      TSPseudoSetMaxTimeStep_Pseudo);
717:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",
718:                     "TSPseudoIncrementDtFromInitialDt_Pseudo",
719:                      TSPseudoIncrementDtFromInitialDt_Pseudo);
720:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C",
721:                     "TSPseudoSetTimeStep_Pseudo",
722:                      TSPseudoSetTimeStep_Pseudo);
723:   return(0);
724: }
725: EXTERN_C_END

729: /*@C
730:    TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping.
731:    Use with TSPseudoSetTimeStep().

733:    Collective on TS

735:    Input Parameters:
736: .  ts - the timestep context
737: .  dtctx - unused timestep context

739:    Output Parameter:
740: .  newdt - the timestep to use for the next step

742:    Level: advanced

744: .keywords: timestep, pseudo, default

746: .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep()
747: @*/
748: PetscErrorCode  TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx)
749: {
750:   TS_Pseudo      *pseudo = (TS_Pseudo*)ts->data;
751:   PetscReal      inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous;

755:   VecZeroEntries(pseudo->xdot);
756:   TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);
757:   VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);
758:   if (pseudo->fnorm_initial == 0.0) {
759:     /* first time through so compute initial function norm */
760:     pseudo->fnorm_initial = pseudo->fnorm;
761:     fnorm_previous        = pseudo->fnorm;
762:   }
763:   if (pseudo->fnorm == 0.0) {
764:     *newdt = 1.e12*inc*ts->time_step;
765:   } else if (pseudo->increment_dt_from_initial_dt) {
766:     *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm;
767:   } else {
768:     *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm;
769:   }
770:   if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max);
771:   pseudo->fnorm_previous = pseudo->fnorm;
772:   return(0);
773: }