Actual source code: beuler.c

  1: #define PETSCTS_DLL

  3: /*
  4:        Code for Timestepping with implicit backwards Euler.
  5: */
 6:  #include private/tsimpl.h

  8: typedef struct {
  9:   Vec  update;      /* work vector where new solution is formed */
 10:   Vec  func;        /* work vector where F(t[i],u[i]) is stored */
 11:   Vec  rhs;         /* work vector for RHS; vec_sol/dt */
 12: } TS_BEuler;

 14: /*------------------------------------------------------------------------------*/
 15: /* 
 16:    Set ts->A = ts->Arhs = 1/dt*Alhs - Arhs, used in KSPSolve() 
 17: */
 20: PetscErrorCode TSSetKSPOperators_BEuler(TS ts)
 21: {
 23:   PetscScalar    mdt = 1.0/ts->time_step;

 26:   if (!ts->A){
 27:     PetscObjectReference((PetscObject)ts->Arhs);
 28:     ts->A = ts->Arhs;
 29:   }

 31:   MatScale(ts->A,-1.0);
 32:   if (ts->Alhs){
 33:     MatAXPY(ts->A,mdt,ts->Alhs,DIFFERENT_NONZERO_PATTERN);
 34:   } else {
 35:     MatShift(ts->A,mdt);
 36:   }
 37:   return(0);
 38: }

 40: /*
 41:     Version for linear PDE where RHS does not depend on time. Has built a
 42:   single matrix that is to be used for all timesteps.
 43: */
 46: static PetscErrorCode TSStep_BEuler_Linear_Constant_Matrix(TS ts,PetscInt *steps,PetscReal *ptime)
 47: {
 48:   TS_BEuler      *beuler = (TS_BEuler*)ts->data;
 49:   Vec            sol = ts->vec_sol,update = beuler->update;
 50:   Vec            rhs = beuler->rhs;
 52:   PetscInt       i,max_steps = ts->max_steps,its;
 53:   PetscScalar    mdt = 1.0/ts->time_step;
 54:   KSP            ksp;

 57:   TSGetKSP(ts,&ksp);
 58:   *steps = -ts->steps;
 59:   TSMonitor(ts,ts->steps,ts->ptime,sol);

 61:   /* set initial guess to be previous solution */
 62:   VecCopy(sol,update);

 64:   for (i=0; i<max_steps; i++) {
 65:     if (ts->ptime + ts->time_step > ts->max_time) break;

 67:     /* set rhs = 1/dt*Alhs*sol */
 68:     if (ts->Alhs){
 69:       MatMult(ts->Alhs,sol,rhs);
 70:     } else {
 71:       VecCopy(sol,rhs);
 72:     }
 73:     VecScale(rhs,mdt);

 75:     ts->ptime += ts->time_step;

 77:     /* solve (1/dt*Alhs - A)*update = rhs */
 78:     KSPSolve(ts->ksp,rhs,update);
 79:     KSPGetIterationNumber(ksp,&its);
 80:     ts->linear_its += its;
 81:     VecCopy(update,sol);
 82:     ts->steps++;
 83:     TSMonitor(ts,ts->steps,ts->ptime,sol);
 84:   }

 86:   *steps += ts->steps;
 87:   *ptime  = ts->ptime;
 88:   return(0);
 89: }

 91: /*
 92:       Version where matrix depends on time 
 93: */
 96: static PetscErrorCode TSStep_BEuler_Linear_Variable_Matrix(TS ts,PetscInt *steps,PetscReal *ptime)
 97: {
 98:   TS_BEuler      *beuler = (TS_BEuler*)ts->data;
 99:   Vec            sol = ts->vec_sol,update = beuler->update,rhs = beuler->rhs;
101:   PetscInt       i,max_steps = ts->max_steps,its;
102:   PetscScalar    mdt = 1.0/ts->time_step;
103:   PetscReal      t_mid;
104:   MatStructure   str;
105:   KSP            ksp;

108:   TSGetKSP(ts,&ksp);
109:   *steps = -ts->steps;
110:   TSMonitor(ts,ts->steps,ts->ptime,sol);

112:   /* set initial guess to be previous solution */
113:   VecCopy(sol,update);

115:   for (i=0; i<max_steps; i++) {
116:     if (ts->ptime +ts->time_step > ts->max_time) break;

118:     /* set rhs = 1/dt*Alhs(t_mid)*sol */
119:     if (ts->Alhs){
120:       /* evaluate lhs matrix function at t_mid */
121:       t_mid = ts->ptime+ts->time_step/2.0;
122:       (*ts->ops->lhsmatrix)(ts,t_mid,&ts->Alhs,PETSC_NULL,&str,ts->jacPlhs);
123:       MatMult(ts->Alhs,sol,rhs);
124:     } else {
125:       VecCopy(sol,rhs);
126:     }
127:     VecScale(rhs,mdt);

129:     ts->ptime += ts->time_step;

131:     /* evaluate rhs matrix function at current ptime */
132:     (*ts->ops->rhsmatrix)(ts,ts->ptime,&ts->Arhs,&ts->B,&str,ts->jacP);

134:     /* set ts->A = ts->Arhs = 1/dt*Alhs - Arhs, used in KSPSolve() */
135:     TSSetKSPOperators_BEuler(ts);
136:     KSPSetOperators(ts->ksp,ts->A,ts->A,SAME_NONZERO_PATTERN);

138:     /* solve (1/dt*Alhs(t_mid) - A(t_n+1))*update = rhs */
139:     KSPSolve(ts->ksp,rhs,update);
140:     KSPGetIterationNumber(ksp,&its);
141:     ts->linear_its += its;
142:     VecCopy(update,sol);
143:     ts->steps++;
144:     TSMonitor(ts,ts->steps,ts->ptime,sol);
145:   }

147:   *steps += ts->steps;
148:   *ptime  = ts->ptime;
149:   return(0);
150: }
151: /*
152:     Version for nonlinear PDE.
153: */
156: static PetscErrorCode TSStep_BEuler_Nonlinear(TS ts,PetscInt *steps,PetscReal *ptime)
157: {
158:   Vec            sol = ts->vec_sol;
160:   PetscInt       i,max_steps = ts->max_steps,its,lits;
161:   TS_BEuler      *beuler = (TS_BEuler*)ts->data;
162: 
164:   *steps = -ts->steps;
165:   TSMonitor(ts,ts->steps,ts->ptime,sol);

167:   for (i=0; i<max_steps; i++) {
168:     if (ts->ptime + ts->time_step > ts->max_time) break;
169:     ts->ptime += ts->time_step;
170:     VecCopy(sol,beuler->update);
171:     SNESSolve(ts->snes,PETSC_NULL,beuler->update);
172:     SNESGetLinearSolveIterations(ts->snes,&lits);
173:     SNESGetIterationNumber(ts->snes,&its);
174:     ts->nonlinear_its += its; ts->linear_its += lits;
175:     VecCopy(beuler->update,sol);
176:     ts->steps++;
177:     TSMonitor(ts,ts->steps,ts->ptime,sol);
178:   }

180:   *steps += ts->steps;
181:   *ptime  = ts->ptime;
182:   return(0);
183: }

185: /*------------------------------------------------------------*/
188: static PetscErrorCode TSDestroy_BEuler(TS ts)
189: {
190:   TS_BEuler      *beuler = (TS_BEuler*)ts->data;

194:   if (beuler->update) {VecDestroy(beuler->update);}
195:   if (beuler->func) {VecDestroy(beuler->func);}
196:   if (beuler->rhs) {VecDestroy(beuler->rhs);}
197:   PetscFree(beuler);
198:   return(0);
199: }

201: /* 
202:     This defines the nonlinear equation that is to be solved with SNES
203:       1/dt* (U^{n+1} - U^{n}) - F(U^{n+1}) 
204: */
207: PetscErrorCode TSBEulerFunction(SNES snes,Vec x,Vec y,void *ctx)
208: {
209:   TS             ts = (TS) ctx;
210:   PetscScalar    mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
212:   PetscInt       i,n;

215:   /* apply user-provided function */
216:   TSComputeRHSFunction(ts,ts->ptime,x,y);
217:   VecGetArray(ts->vec_sol,&un);
218:   VecGetArray(x,&unp1);
219:   VecGetArray(y,&Funp1);
220:   VecGetLocalSize(x,&n);
221:   for (i=0; i<n; i++) {
222:     Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
223:   }
224:   VecRestoreArray(ts->vec_sol,&un);
225:   VecRestoreArray(x,&unp1);
226:   VecRestoreArray(y,&Funp1);
227:   return(0);
228: }

230: /*
231:    This constructs the Jacobian needed for SNES 
232:      J = I/dt - J_{F}   where J_{F} is the given Jacobian of F at t_{n+1}.
233:      x  - input vector
234:      AA - Jacobian matrix 
235:      BB - preconditioner matrix, usually the same as AA
236: */
239: PetscErrorCode TSBEulerJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
240: {
241:   TS             ts = (TS) ctx;

245:   /* construct user's Jacobian */
246:   TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);

248:   /* shift and scale Jacobian */
249:   /* this test is a undesirable hack, we assume that if it is MATMFFD then it is
250:      obtained from -snes_mf_operator and there is computed directly from the 
251:      FormFunction() SNES is given and therefor does not need to be shifted/scaled
252:      BUT maybe it could be MATMFFD and does require shift in some other case??? */
253:   TSSetKSPOperators_BEuler(ts);
254:   return(0);
255: }

257: /* ------------------------------------------------------------*/
260: static PetscErrorCode TSSetUp_BEuler_Linear_Constant_Matrix(TS ts)
261: {
262:   TS_BEuler      *beuler = (TS_BEuler*)ts->data;

266:   VecDuplicate(ts->vec_sol,&beuler->update);
267:   VecDuplicate(ts->vec_sol,&beuler->rhs);
268: 
269:   /* build linear system to be solved - should move into TSStep() if dt changes! */
270:   /* Set ts->A = ts->Arhs = 1/dt*Alhs - Arhs, used in KSPSolve() */
271:   TSSetKSPOperators_BEuler(ts);
272:   KSPSetOperators(ts->ksp,ts->A,ts->A,SAME_NONZERO_PATTERN);
273:   return(0);
274: }

278: static PetscErrorCode TSSetUp_BEuler_Linear_Variable_Matrix(TS ts)
279: {
280:   TS_BEuler      *beuler = (TS_BEuler*)ts->data;

284:   VecDuplicate(ts->vec_sol,&beuler->update);
285:   VecDuplicate(ts->vec_sol,&beuler->rhs);
286:   return(0);
287: }

291: static PetscErrorCode TSSetUp_BEuler_Nonlinear(TS ts)
292: {
293:   TS_BEuler      *beuler = (TS_BEuler*)ts->data;

297:   VecDuplicate(ts->vec_sol,&beuler->update);
298:   VecDuplicate(ts->vec_sol,&beuler->func);
299:   SNESSetFunction(ts->snes,beuler->func,TSBEulerFunction,ts);
300:   SNESSetJacobian(ts->snes,ts->Arhs,ts->B,TSBEulerJacobian,ts);
301:   return(0);
302: }
303: /*------------------------------------------------------------*/

307: static PetscErrorCode TSSetFromOptions_BEuler_Linear(TS ts)
308: {
310:   return(0);
311: }

315: static PetscErrorCode TSSetFromOptions_BEuler_Nonlinear(TS ts)
316: {
318:   return(0);
319: }

323: static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
324: {
326:   return(0);
327: }

329: /* ------------------------------------------------------------ */
330: /*MC
331:       TS_BEULER - ODE solver using the implicit backward Euler method

333:   Level: beginner

335: .seealso:  TSCreate(), TS, TSSetType(), TS_EULER

337: M*/
341: PetscErrorCode  TSCreate_BEuler(TS ts)
342: {
343:   TS_BEuler      *beuler;

347:   ts->ops->destroy = TSDestroy_BEuler;
348:   ts->ops->view    = TSView_BEuler;

350:   if (ts->problem_type == TS_LINEAR) {
351:     if (!ts->Arhs) {
352:       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set rhs matrix for linear problem");
353:     }
354:     if (!ts->ops->rhsmatrix) {
355:       ts->ops->setup  = TSSetUp_BEuler_Linear_Constant_Matrix;
356:       ts->ops->step   = TSStep_BEuler_Linear_Constant_Matrix;
357:     } else {
358:       ts->ops->setup  = TSSetUp_BEuler_Linear_Variable_Matrix;
359:       ts->ops->step   = TSStep_BEuler_Linear_Variable_Matrix;
360:     }
361:     ts->ops->setfromoptions  = TSSetFromOptions_BEuler_Linear;
362:     KSPCreate(((PetscObject)ts)->comm,&ts->ksp);
363:     PetscObjectIncrementTabLevel((PetscObject)ts->ksp,(PetscObject)ts,1);
364:     KSPSetInitialGuessNonzero(ts->ksp,PETSC_TRUE);
365:   } else if (ts->problem_type == TS_NONLINEAR) {
366:     ts->ops->setup           = TSSetUp_BEuler_Nonlinear;
367:     ts->ops->step            = TSStep_BEuler_Nonlinear;
368:     ts->ops->setfromoptions  = TSSetFromOptions_BEuler_Nonlinear;
369:     SNESCreate(((PetscObject)ts)->comm,&ts->snes);
370:     PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);
371:   } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"No such problem");

373:   PetscNewLog(ts,TS_BEuler,&beuler);
374:   ts->data = (void*)beuler;

376:   return(0);
377: }