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: */
19: static int TSStep_BEuler_Linear_Constant_Matrix(TS ts,int *steps,PetscReal *ptime)
20: {
21: TS_BEuler *beuler = (TS_BEuler*)ts->data;
22: Vec sol = ts->vec_sol,update = beuler->update;
23: Vec rhs = beuler->rhs;
24: int ierr,i,max_steps = ts->max_steps,its;
25: PetscScalar mdt = 1.0/ts->time_step;
26:
28: *steps = -ts->steps;
29: TSMonitor(ts,ts->steps,ts->ptime,sol);
31: /* set initial guess to be previous solution */
32: VecCopy(sol,update);
34: for (i=0; i<max_steps; i++) {
35: VecCopy(sol,rhs);
36: VecScale(&mdt,rhs);
37: /* apply user-provided boundary conditions (only needed if they are time dependent) */
38: TSComputeRHSBoundaryConditions(ts,ts->ptime,rhs);
40: ts->ptime += ts->time_step;
41: if (ts->ptime > ts->max_time) break;
42: SLESSolve(ts->sles,rhs,update,&its);
43: ts->linear_its += its;
44: VecCopy(update,sol);
45: ts->steps++;
46: TSMonitor(ts,ts->steps,ts->ptime,sol);
47: }
49: *steps += ts->steps;
50: *ptime = ts->ptime;
51: return(0);
52: }
53: /*
54: Version where matrix depends on time
55: */
56: static int TSStep_BEuler_Linear_Variable_Matrix(TS ts,int *steps,PetscReal *ptime)
57: {
58: TS_BEuler *beuler = (TS_BEuler*)ts->data;
59: Vec sol = ts->vec_sol,update = beuler->update,rhs = beuler->rhs;
60: int ierr,i,max_steps = ts->max_steps,its;
61: PetscScalar mdt = 1.0/ts->time_step,mone = -1.0;
62: MatStructure str;
65: *steps = -ts->steps;
66: TSMonitor(ts,ts->steps,ts->ptime,sol);
68: /* set initial guess to be previous solution */
69: VecCopy(sol,update);
71: for (i=0; i<max_steps; i++) {
72: VecCopy(sol,rhs);
73: VecScale(&mdt,rhs);
74: /* apply user-provided boundary conditions (only needed if they are time dependent) */
75: TSComputeRHSBoundaryConditions(ts,ts->ptime,rhs);
77: ts->ptime += ts->time_step;
78: if (ts->ptime > ts->max_time) break;
79: /*
80: evaluate matrix function
81: */
82: (*ts->ops->rhsmatrix)(ts,ts->ptime,&ts->A,&ts->B,&str,ts->jacP);
83: MatScale(&mone,ts->A);
84: MatShift(&mdt,ts->A);
85: if (ts->B != ts->A && str != SAME_PRECONDITIONER) {
86: MatScale(&mone,ts->B);
87: MatShift(&mdt,ts->B);
88: }
89: SLESSetOperators(ts->sles,ts->A,ts->B,str);
90: SLESSolve(ts->sles,rhs,update,&its);
91: ts->linear_its += its;
92: VecCopy(update,sol);
93: ts->steps++;
94: TSMonitor(ts,ts->steps,ts->ptime,sol);
95: }
97: *steps += ts->steps;
98: *ptime = ts->ptime;
99: return(0);
100: }
101: /*
102: Version for nonlinear PDE.
103: */
104: static int TSStep_BEuler_Nonlinear(TS ts,int *steps,PetscReal *ptime)
105: {
106: Vec sol = ts->vec_sol;
107: int ierr,i,max_steps = ts->max_steps,its,lits;
108: TS_BEuler *beuler = (TS_BEuler*)ts->data;
109:
111: *steps = -ts->steps;
112: TSMonitor(ts,ts->steps,ts->ptime,sol);
114: for (i=0; i<max_steps; i++) {
115: ts->ptime += ts->time_step;
116: if (ts->ptime > ts->max_time) break;
117: VecCopy(sol,beuler->update);
118: SNESSolve(ts->snes,beuler->update,&its);
119: SNESGetNumberLinearIterations(ts->snes,&lits);
120: ts->nonlinear_its += its; ts->linear_its += lits;
121: VecCopy(beuler->update,sol);
122: ts->steps++;
123: TSMonitor(ts,ts->steps,ts->ptime,sol);
124: }
126: *steps += ts->steps;
127: *ptime = ts->ptime;
128: return(0);
129: }
131: /*------------------------------------------------------------*/
132: static int TSDestroy_BEuler(TS ts)
133: {
134: TS_BEuler *beuler = (TS_BEuler*)ts->data;
135: int ierr;
138: if (beuler->update) {VecDestroy(beuler->update);}
139: if (beuler->func) {VecDestroy(beuler->func);}
140: if (beuler->rhs) {VecDestroy(beuler->rhs);}
141: PetscFree(beuler);
142: return(0);
143: }
147: /*
148: This defines the nonlinear equation that is to be solved with SNES
150: U^{n+1} - dt*F(U^{n+1}) - U^{n}
151: */
152: int TSBEulerFunction(SNES snes,Vec x,Vec y,void *ctx)
153: {
154: TS ts = (TS) ctx;
155: PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
156: int ierr,i,n;
159: /* apply user-provided function */
160: TSComputeRHSFunction(ts,ts->ptime,x,y);
161: /* (u^{n+1} - U^{n})/dt - F(u^{n+1}) */
162: VecGetArray(ts->vec_sol,&un);
163: VecGetArray(x,&unp1);
164: VecGetArray(y,&Funp1);
165: VecGetLocalSize(x,&n);
167: for (i=0; i<n; i++) {
168: Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
169: }
170: VecRestoreArray(ts->vec_sol,&un);
171: VecRestoreArray(x,&unp1);
172: VecRestoreArray(y,&Funp1);
173: return(0);
174: }
176: /*
177: This constructs the Jacobian needed for SNES
179: J = I/dt - J_{F} where J_{F} is the given Jacobian of F.
180: */
181: int TSBEulerJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
182: {
183: TS ts = (TS) ctx;
184: int ierr;
185: PetscScalar mone = -1.0,mdt = 1.0/ts->time_step;
188: /* construct user's Jacobian */
189: TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);
191: /* shift and scale Jacobian */
192: MatScale(&mone,*AA);
193: MatShift(&mdt,*AA);
194: if (*BB != *AA && *str != SAME_PRECONDITIONER) {
195: MatScale(&mone,*BB);
196: MatShift(&mdt,*BB);
197: }
199: return(0);
200: }
202: /* ------------------------------------------------------------*/
203: static int TSSetUp_BEuler_Linear_Constant_Matrix(TS ts)
204: {
205: TS_BEuler *beuler = (TS_BEuler*)ts->data;
206: int ierr;
207: PetscScalar mdt = 1.0/ts->time_step,mone = -1.0;
210: SLESSetFromOptions(ts->sles);
211: VecDuplicate(ts->vec_sol,&beuler->update);
212: VecDuplicate(ts->vec_sol,&beuler->rhs);
213:
214: /* build linear system to be solved */
215: MatScale(&mone,ts->A);
216: MatShift(&mdt,ts->A);
217: if (ts->A != ts->B) {
218: MatScale(&mone,ts->B);
219: MatShift(&mdt,ts->B);
220: }
221: SLESSetOperators(ts->sles,ts->A,ts->B,SAME_NONZERO_PATTERN);
222: return(0);
223: }
225: static int TSSetUp_BEuler_Linear_Variable_Matrix(TS ts)
226: {
227: TS_BEuler *beuler = (TS_BEuler*)ts->data;
228: int ierr;
231: SLESSetFromOptions(ts->sles);
232: VecDuplicate(ts->vec_sol,&beuler->update);
233: VecDuplicate(ts->vec_sol,&beuler->rhs);
234: return(0);
235: }
237: static int TSSetUp_BEuler_Nonlinear(TS ts)
238: {
239: TS_BEuler *beuler = (TS_BEuler*)ts->data;
240: int ierr;
243: VecDuplicate(ts->vec_sol,&beuler->update);
244: VecDuplicate(ts->vec_sol,&beuler->func);
245: SNESSetFunction(ts->snes,beuler->func,TSBEulerFunction,ts);
246: SNESSetJacobian(ts->snes,ts->A,ts->B,TSBEulerJacobian,ts);
247: return(0);
248: }
249: /*------------------------------------------------------------*/
251: static int TSSetFromOptions_BEuler_Linear(TS ts)
252: {
254:
255: return(0);
256: }
258: static int TSSetFromOptions_BEuler_Nonlinear(TS ts)
259: {
261:
262: return(0);
263: }
265: static int TSView_BEuler(TS ts,PetscViewer viewer)
266: {
268: return(0);
269: }
271: /* ------------------------------------------------------------ */
272: EXTERN_C_BEGIN
273: int TSCreate_BEuler(TS ts)
274: {
275: TS_BEuler *beuler;
276: int ierr;
277: KSP ksp;
280: ts->ops->destroy = TSDestroy_BEuler;
281: ts->ops->view = TSView_BEuler;
283: if (ts->problem_type == TS_LINEAR) {
284: if (!ts->A) {
285: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set rhs matrix for linear problem");
286: }
287: if (!ts->ops->rhsmatrix) {
288: ts->ops->setup = TSSetUp_BEuler_Linear_Constant_Matrix;
289: ts->ops->step = TSStep_BEuler_Linear_Constant_Matrix;
290: } else {
291: ts->ops->setup = TSSetUp_BEuler_Linear_Variable_Matrix;
292: ts->ops->step = TSStep_BEuler_Linear_Variable_Matrix;
293: }
294: ts->ops->setfromoptions = TSSetFromOptions_BEuler_Linear;
295: SLESCreate(ts->comm,&ts->sles);
296: SLESGetKSP(ts->sles,&ksp);
297: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
298: } else if (ts->problem_type == TS_NONLINEAR) {
299: ts->ops->setup = TSSetUp_BEuler_Nonlinear;
300: ts->ops->step = TSStep_BEuler_Nonlinear;
301: ts->ops->setfromoptions = TSSetFromOptions_BEuler_Nonlinear;
302: SNESCreate(ts->comm,&ts->snes);
303: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"No such problem");
305: PetscNew(TS_BEuler,&beuler);
306: PetscLogObjectMemory(ts,sizeof(TS_BEuler));
307: ierr = PetscMemzero(beuler,sizeof(TS_BEuler));
308: ts->data = (void*)beuler;
310: return(0);
311: }
312: EXTERN_C_END