Actual source code: theta.c

  1: #define PETSCTS_DLL

  3: /*
  4:   Code for timestepping with implicit Theta method

  6:   Notes:
  7:   This method can be applied to DAE.

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

 11:   Theta | Theta
 12:   -------------
 13:         |  1

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

 17:   X_i = x + h sum_j a_ij X'_j

 19:   is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i)
 20: */
 21:  #include private/tsimpl.h

 23: typedef struct {
 24:   Vec X,Xdot;                   /* Storage for one stage */
 25:   Vec res;                      /* DAE residuals */
 26:   PetscTruth extrapolate;
 27:   PetscReal Theta;
 28:   PetscReal shift;
 29:   PetscReal stage_time;
 30: } TS_Theta;

 34: static PetscErrorCode TSStep_Theta(TS ts,PetscInt *steps,PetscReal *ptime)
 35: {
 36:   TS_Theta       *th = (TS_Theta*)ts->data;
 37:   PetscInt       i,max_steps = ts->max_steps,its,lits;

 41:   *steps = -ts->steps;
 42:   *ptime = ts->ptime;

 44:   TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);

 46:   for (i=0; i<max_steps; i++) {
 47:     if (ts->ptime + ts->time_step > ts->max_time) break;
 48:     TSPreStep(ts);

 50:     th->stage_time = ts->ptime + th->Theta*ts->time_step;
 51:     th->shift = 1./(th->Theta*ts->time_step);

 53:     if (th->extrapolate) {
 54:       VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);
 55:     } else {
 56:       VecCopy(ts->vec_sol,th->X);
 57:     }
 58:     SNESSolve(ts->snes,PETSC_NULL,th->X);
 59:     SNESGetIterationNumber(ts->snes,&its);
 60:     SNESGetLinearSolveIterations(ts->snes,&lits);
 61:     ts->nonlinear_its += its; ts->linear_its += lits;

 63:     VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);
 64:     VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);
 65:     ts->ptime += ts->time_step;
 66:     ts->steps++;

 68:     TSPostStep(ts);
 69:     TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
 70:   }

 72:   *steps += ts->steps;
 73:   *ptime  = ts->ptime;
 74:   return(0);
 75: }

 77: /*------------------------------------------------------------*/
 80: static PetscErrorCode TSDestroy_Theta(TS ts)
 81: {
 82:   TS_Theta       *th = (TS_Theta*)ts->data;
 83:   PetscErrorCode  ierr;

 86:   if (th->X)    {VecDestroy(th->X);}
 87:   if (th->Xdot) {VecDestroy(th->Xdot);}
 88:   if (th->res)  {VecDestroy(th->res);}
 89:   PetscFree(th);
 90:   return(0);
 91: }

 93: /*
 94:   This defines the nonlinear equation that is to be solved with SNES
 95:   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
 96: */
 99: static PetscErrorCode TSThetaFunction(SNES snes,Vec x,Vec y,void *ctx)
100: {
101:   TS        ts = (TS)ctx;
102:   TS_Theta *th = (TS_Theta*)ts->data;

106:   VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);
107:   TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y);
108:   return(0);
109: }

113: static PetscErrorCode TSThetaJacobian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,void *ctx)
114: {
115:   TS        ts = (TS)ctx;
116:   TS_Theta *th = (TS_Theta*)ts->data;

120:   /* th->Xdot has already been computed in TSThetaFunction (SNES guarantees this) */
121:   TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str);
122:   return(0);
123: }


128: static PetscErrorCode TSSetUp_Theta(TS ts)
129: {
130:   TS_Theta *th = (TS_Theta*)ts->data;

134:   if (ts->problem_type == TS_LINEAR) {
135:     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
136:   }
137:   VecDuplicate(ts->vec_sol,&th->X);
138:   VecDuplicate(ts->vec_sol,&th->Xdot);
139:   VecDuplicate(ts->vec_sol,&th->res);
140:   SNESSetFunction(ts->snes,th->res,TSThetaFunction,ts);
141:   /* This is nasty.  SNESSetFromOptions() is usually called in TSSetFromOptions().  With -snes_mf_operator, it will
142:   replace A and we don't want to mess with that.  With -snes_mf, A and B will be replaced as well as the function and
143:   context.  Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets
144:   the Jacobian user context to snes->funP, it will actually be NULL.  This is not a problem because both snes->funP and
145:   snes->jacP should be the TS. */
146:   {
147:     Mat A,B;
148:     PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
149:     void *ctx;
150:     SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);
151:     SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&TSThetaJacobian,ctx?ctx:ts);
152:   }
153:   return(0);
154: }
155: /*------------------------------------------------------------*/

159: static PetscErrorCode TSSetFromOptions_Theta(TS ts)
160: {
161:   TS_Theta *th = (TS_Theta*)ts->data;

165:   PetscOptionsHead("Theta ODE solver options");
166:   {
167:     PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);
168:     PetscOptionsTruth("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);
169:   }
170:   PetscOptionsTail();
171:   return(0);
172: }

176: static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
177: {
178:   TS_Theta       *th = (TS_Theta*)ts->data;
179:   PetscTruth      iascii;
180:   PetscErrorCode  ierr;

183:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
184:   if (iascii) {
185:     PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);
186:     PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");
187:   } else {
188:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name);
189:   }
190:   return(0);
191: }

193: /* ------------------------------------------------------------ */
194: /*MC
195:       TSTHETA - DAE solver using the implicit Theta method

197:   Level: beginner

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

201: M*/
205: PetscErrorCode  TSCreate_Theta(TS ts)
206: {
207:   TS_Theta       *th;

211:   ts->ops->destroy        = TSDestroy_Theta;
212:   ts->ops->view           = TSView_Theta;
213:   ts->ops->setup          = TSSetUp_Theta;
214:   ts->ops->step           = TSStep_Theta;
215:   ts->ops->setfromoptions = TSSetFromOptions_Theta;

217:   ts->problem_type = TS_NONLINEAR;
218:   SNESCreate(((PetscObject)ts)->comm,&ts->snes);
219:   PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);

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

224:   th->extrapolate = PETSC_TRUE;
225:   th->Theta       = 0.5;

227:   return(0);
228: }