Actual source code: beuler.c
1: /*$Id: beuler.c,v 1.61 2001/09/11 16:34:19 bsmith Exp $*/
2: /*
3: Code for Timestepping with implicit backwards Euler.
4: */
5: #include src/ts/tsimpl.h
7: typedef struct {
8: Vec update; /* work vector where new solution is formed */
9: Vec func; /* work vector where F(t[i],u[i]) is stored */
10: Vec rhs; /* work vector for RHS; vec_sol/dt */
11: } TS_BEuler;
13: /*------------------------------------------------------------------------------*/
15: /*
16: Version for linear PDE where RHS does not depend on time. Has built a
17: single matrix that is to be used for all timesteps.
18: */
21: static int TSStep_BEuler_Linear_Constant_Matrix(TS ts,int *steps,PetscReal *ptime)
22: {
23: TS_BEuler *beuler = (TS_BEuler*)ts->data;
24: Vec sol = ts->vec_sol,update = beuler->update;
25: Vec rhs = beuler->rhs;
26: int ierr,i,max_steps = ts->max_steps,its;
27: PetscScalar mdt = 1.0/ts->time_step;
28: KSP ksp;
31: TSGetKSP(ts,&ksp);
32: *steps = -ts->steps;
33: TSMonitor(ts,ts->steps,ts->ptime,sol);
35: /* set initial guess to be previous solution */
36: VecCopy(sol,update);
38: for (i=0; i<max_steps; i++) {
39: VecCopy(sol,rhs);
40: VecScale(&mdt,rhs);
41: /* apply user-provided boundary conditions (only needed if they are time dependent) */
42: TSComputeRHSBoundaryConditions(ts,ts->ptime,rhs);
44: ts->ptime += ts->time_step;
45: if (ts->ptime > ts->max_time) break;
46: KSPSetRhs(ts->ksp,rhs);
47: KSPSetSolution(ts->ksp,update);
48: KSPSolve(ts->ksp);
49: KSPGetIterationNumber(ksp,&its);
50: ts->linear_its += its;
51: VecCopy(update,sol);
52: ts->steps++;
53: TSMonitor(ts,ts->steps,ts->ptime,sol);
54: }
56: *steps += ts->steps;
57: *ptime = ts->ptime;
58: return(0);
59: }
60: /*
61: Version where matrix depends on time
62: */
65: static int TSStep_BEuler_Linear_Variable_Matrix(TS ts,int *steps,PetscReal *ptime)
66: {
67: TS_BEuler *beuler = (TS_BEuler*)ts->data;
68: Vec sol = ts->vec_sol,update = beuler->update,rhs = beuler->rhs;
69: int ierr,i,max_steps = ts->max_steps,its;
70: PetscScalar mdt = 1.0/ts->time_step,mone = -1.0;
71: MatStructure str;
72: KSP ksp;
75: TSGetKSP(ts,&ksp);
76: *steps = -ts->steps;
77: TSMonitor(ts,ts->steps,ts->ptime,sol);
79: /* set initial guess to be previous solution */
80: VecCopy(sol,update);
82: for (i=0; i<max_steps; i++) {
83: VecCopy(sol,rhs);
84: VecScale(&mdt,rhs);
85: /* apply user-provided boundary conditions (only needed if they are time dependent) */
86: TSComputeRHSBoundaryConditions(ts,ts->ptime,rhs);
88: ts->ptime += ts->time_step;
89: if (ts->ptime > ts->max_time) break;
90: /*
91: evaluate matrix function
92: */
93: (*ts->ops->rhsmatrix)(ts,ts->ptime,&ts->A,&ts->B,&str,ts->jacP);
94: MatScale(&mone,ts->A);
95: MatShift(&mdt,ts->A);
96: if (ts->B != ts->A && str != SAME_PRECONDITIONER) {
97: MatScale(&mone,ts->B);
98: MatShift(&mdt,ts->B);
99: }
100: KSPSetOperators(ts->ksp,ts->A,ts->B,str);
101: KSPSetRhs(ts->ksp,rhs);
102: KSPSetSolution(ts->ksp,update);
103: KSPSolve(ts->ksp);
104: KSPGetIterationNumber(ksp,&its);
105: ts->linear_its += its;
106: VecCopy(update,sol);
107: ts->steps++;
108: TSMonitor(ts,ts->steps,ts->ptime,sol);
109: }
111: *steps += ts->steps;
112: *ptime = ts->ptime;
113: return(0);
114: }
115: /*
116: Version for nonlinear PDE.
117: */
120: static int TSStep_BEuler_Nonlinear(TS ts,int *steps,PetscReal *ptime)
121: {
122: Vec sol = ts->vec_sol;
123: int ierr,i,max_steps = ts->max_steps,its,lits;
124: TS_BEuler *beuler = (TS_BEuler*)ts->data;
125:
127: *steps = -ts->steps;
128: TSMonitor(ts,ts->steps,ts->ptime,sol);
130: for (i=0; i<max_steps; i++) {
131: ts->ptime += ts->time_step;
132: if (ts->ptime > ts->max_time) break;
133: VecCopy(sol,beuler->update);
134: SNESSolve(ts->snes,beuler->update);
135: SNESGetNumberLinearIterations(ts->snes,&lits);
136: SNESGetIterationNumber(ts->snes,&its);
137: ts->nonlinear_its += its; ts->linear_its += lits;
138: VecCopy(beuler->update,sol);
139: ts->steps++;
140: TSMonitor(ts,ts->steps,ts->ptime,sol);
141: }
143: *steps += ts->steps;
144: *ptime = ts->ptime;
145: return(0);
146: }
148: /*------------------------------------------------------------*/
151: static int TSDestroy_BEuler(TS ts)
152: {
153: TS_BEuler *beuler = (TS_BEuler*)ts->data;
154: int ierr;
157: if (beuler->update) {VecDestroy(beuler->update);}
158: if (beuler->func) {VecDestroy(beuler->func);}
159: if (beuler->rhs) {VecDestroy(beuler->rhs);}
160: PetscFree(beuler);
161: return(0);
162: }
166: /*
167: This defines the nonlinear equation that is to be solved with SNES
169: U^{n+1} - dt*F(U^{n+1}) - U^{n}
170: */
173: int TSBEulerFunction(SNES snes,Vec x,Vec y,void *ctx)
174: {
175: TS ts = (TS) ctx;
176: PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
177: int ierr,i,n;
180: /* apply user-provided function */
181: TSComputeRHSFunction(ts,ts->ptime,x,y);
182: /* (u^{n+1} - U^{n})/dt - F(u^{n+1}) */
183: VecGetArray(ts->vec_sol,&un);
184: VecGetArray(x,&unp1);
185: VecGetArray(y,&Funp1);
186: VecGetLocalSize(x,&n);
188: for (i=0; i<n; i++) {
189: Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
190: }
191: VecRestoreArray(ts->vec_sol,&un);
192: VecRestoreArray(x,&unp1);
193: VecRestoreArray(y,&Funp1);
194: return(0);
195: }
197: /*
198: This constructs the Jacobian needed for SNES
200: J = I/dt - J_{F} where J_{F} is the given Jacobian of F.
201: */
204: int TSBEulerJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
205: {
206: TS ts = (TS) ctx;
207: int ierr;
208: PetscScalar mone = -1.0,mdt = 1.0/ts->time_step;
211: /* construct user's Jacobian */
212: TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);
214: /* shift and scale Jacobian */
215: MatScale(&mone,*AA);
216: MatShift(&mdt,*AA);
217: if (*BB != *AA && *str != SAME_PRECONDITIONER) {
218: MatScale(&mone,*BB);
219: MatShift(&mdt,*BB);
220: }
222: return(0);
223: }
225: /* ------------------------------------------------------------*/
228: static int TSSetUp_BEuler_Linear_Constant_Matrix(TS ts)
229: {
230: TS_BEuler *beuler = (TS_BEuler*)ts->data;
231: int ierr;
232: PetscScalar mdt = 1.0/ts->time_step,mone = -1.0;
235: KSPSetFromOptions(ts->ksp);
236: VecDuplicate(ts->vec_sol,&beuler->update);
237: VecDuplicate(ts->vec_sol,&beuler->rhs);
238:
239: /* build linear system to be solved */
240: MatScale(&mone,ts->A);
241: MatShift(&mdt,ts->A);
242: if (ts->A != ts->B) {
243: MatScale(&mone,ts->B);
244: MatShift(&mdt,ts->B);
245: }
246: KSPSetOperators(ts->ksp,ts->A,ts->B,SAME_NONZERO_PATTERN);
247: return(0);
248: }
252: static int TSSetUp_BEuler_Linear_Variable_Matrix(TS ts)
253: {
254: TS_BEuler *beuler = (TS_BEuler*)ts->data;
255: int ierr;
258: KSPSetFromOptions(ts->ksp);
259: VecDuplicate(ts->vec_sol,&beuler->update);
260: VecDuplicate(ts->vec_sol,&beuler->rhs);
261: return(0);
262: }
266: static int TSSetUp_BEuler_Nonlinear(TS ts)
267: {
268: TS_BEuler *beuler = (TS_BEuler*)ts->data;
269: int ierr;
272: VecDuplicate(ts->vec_sol,&beuler->update);
273: VecDuplicate(ts->vec_sol,&beuler->func);
274: SNESSetFunction(ts->snes,beuler->func,TSBEulerFunction,ts);
275: SNESSetJacobian(ts->snes,ts->A,ts->B,TSBEulerJacobian,ts);
276: return(0);
277: }
278: /*------------------------------------------------------------*/
282: static int TSSetFromOptions_BEuler_Linear(TS ts)
283: {
285:
286: return(0);
287: }
291: static int TSSetFromOptions_BEuler_Nonlinear(TS ts)
292: {
294:
295: return(0);
296: }
300: static int TSView_BEuler(TS ts,PetscViewer viewer)
301: {
303: return(0);
304: }
306: /* ------------------------------------------------------------ */
307: /*MC
308: TS_BEULER - ODE solver using the implicit backward Euler method
310: .seealso: TSCreate(), TS, TSSetType(), TS_EULER
312: M*/
313: EXTERN_C_BEGIN
316: int TSCreate_BEuler(TS ts)
317: {
318: TS_BEuler *beuler;
319: int ierr;
322: ts->ops->destroy = TSDestroy_BEuler;
323: ts->ops->view = TSView_BEuler;
325: if (ts->problem_type == TS_LINEAR) {
326: if (!ts->A) {
327: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set rhs matrix for linear problem");
328: }
329: if (!ts->ops->rhsmatrix) {
330: ts->ops->setup = TSSetUp_BEuler_Linear_Constant_Matrix;
331: ts->ops->step = TSStep_BEuler_Linear_Constant_Matrix;
332: } else {
333: ts->ops->setup = TSSetUp_BEuler_Linear_Variable_Matrix;
334: ts->ops->step = TSStep_BEuler_Linear_Variable_Matrix;
335: }
336: ts->ops->setfromoptions = TSSetFromOptions_BEuler_Linear;
337: KSPCreate(ts->comm,&ts->ksp);
338: KSPSetInitialGuessNonzero(ts->ksp,PETSC_TRUE);
339: } else if (ts->problem_type == TS_NONLINEAR) {
340: ts->ops->setup = TSSetUp_BEuler_Nonlinear;
341: ts->ops->step = TSStep_BEuler_Nonlinear;
342: ts->ops->setfromoptions = TSSetFromOptions_BEuler_Nonlinear;
343: SNESCreate(ts->comm,&ts->snes);
344: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"No such problem");
346: PetscNew(TS_BEuler,&beuler);
347: PetscLogObjectMemory(ts,sizeof(TS_BEuler));
348: PetscMemzero(beuler,sizeof(TS_BEuler));
349: ts->data = (void*)beuler;
351: return(0);
352: }
353: EXTERN_C_END