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: TSPreStep(ts);
68: /* set rhs = 1/dt*Alhs*sol */
69: if (ts->Alhs){
70: MatMult(ts->Alhs,sol,rhs);
71: } else {
72: VecCopy(sol,rhs);
73: }
74: VecScale(rhs,mdt);
76: ts->ptime += ts->time_step;
78: /* solve (1/dt*Alhs - A)*update = rhs */
79: KSPSolve(ts->ksp,rhs,update);
80: KSPGetIterationNumber(ksp,&its);
81: ts->linear_its += its;
82: VecCopy(update,sol);
83: ts->steps++;
84: TSPostStep(ts);
85: TSMonitor(ts,ts->steps,ts->ptime,sol);
86: }
88: *steps += ts->steps;
89: *ptime = ts->ptime;
90: return(0);
91: }
93: /*
94: Version where matrix depends on time
95: */
98: static PetscErrorCode TSStep_BEuler_Linear_Variable_Matrix(TS ts,PetscInt *steps,PetscReal *ptime)
99: {
100: TS_BEuler *beuler = (TS_BEuler*)ts->data;
101: Vec sol = ts->vec_sol,update = beuler->update,rhs = beuler->rhs;
103: PetscInt i,max_steps = ts->max_steps,its;
104: PetscScalar mdt = 1.0/ts->time_step;
105: PetscReal t_mid;
106: MatStructure str;
107: KSP ksp;
110: TSGetKSP(ts,&ksp);
111: *steps = -ts->steps;
112: TSMonitor(ts,ts->steps,ts->ptime,sol);
114: /* set initial guess to be previous solution */
115: VecCopy(sol,update);
117: for (i=0; i<max_steps; i++) {
118: if (ts->ptime +ts->time_step > ts->max_time) break;
120: TSPreStep(ts);
121: /* set rhs = 1/dt*Alhs(t_mid)*sol */
122: if (ts->Alhs){
123: /* evaluate lhs matrix function at t_mid */
124: t_mid = ts->ptime+ts->time_step/2.0;
125: (*ts->ops->lhsmatrix)(ts,t_mid,&ts->Alhs,PETSC_NULL,&str,ts->jacPlhs);
126: MatMult(ts->Alhs,sol,rhs);
127: } else {
128: VecCopy(sol,rhs);
129: }
130: VecScale(rhs,mdt);
132: ts->ptime += ts->time_step;
134: /* evaluate rhs matrix function at current ptime */
135: (*ts->ops->rhsmatrix)(ts,ts->ptime,&ts->Arhs,&ts->B,&str,ts->jacP);
137: /* set ts->A = ts->Arhs = 1/dt*Alhs - Arhs, used in KSPSolve() */
138: TSSetKSPOperators_BEuler(ts);
139: KSPSetOperators(ts->ksp,ts->A,ts->A,SAME_NONZERO_PATTERN);
141: /* solve (1/dt*Alhs(t_mid) - A(t_n+1))*update = rhs */
142: KSPSolve(ts->ksp,rhs,update);
143: KSPGetIterationNumber(ksp,&its);
144: ts->linear_its += its;
145: VecCopy(update,sol);
146: ts->steps++;
147: TSPostStep(ts);
148: TSMonitor(ts,ts->steps,ts->ptime,sol);
149: }
151: *steps += ts->steps;
152: *ptime = ts->ptime;
153: return(0);
154: }
155: /*
156: Version for nonlinear PDE.
157: */
160: static PetscErrorCode TSStep_BEuler_Nonlinear(TS ts,PetscInt *steps,PetscReal *ptime)
161: {
162: Vec sol = ts->vec_sol;
164: PetscInt i,max_steps = ts->max_steps,its,lits;
165: TS_BEuler *beuler = (TS_BEuler*)ts->data;
166:
168: *steps = -ts->steps;
169: TSMonitor(ts,ts->steps,ts->ptime,sol);
171: for (i=0; i<max_steps; i++) {
172: if (ts->ptime + ts->time_step > ts->max_time) break;
173: TSPreStep(ts);
174: ts->ptime += ts->time_step;
175: VecCopy(sol,beuler->update);
176: SNESSolve(ts->snes,PETSC_NULL,beuler->update);
177: SNESGetLinearSolveIterations(ts->snes,&lits);
178: SNESGetIterationNumber(ts->snes,&its);
179: ts->nonlinear_its += its; ts->linear_its += lits;
180: VecCopy(beuler->update,sol);
181: ts->steps++;
182: TSPostStep(ts);
183: TSMonitor(ts,ts->steps,ts->ptime,sol);
184: }
186: *steps += ts->steps;
187: *ptime = ts->ptime;
188: return(0);
189: }
191: /*------------------------------------------------------------*/
194: static PetscErrorCode TSDestroy_BEuler(TS ts)
195: {
196: TS_BEuler *beuler = (TS_BEuler*)ts->data;
200: if (beuler->update) {VecDestroy(beuler->update);}
201: if (beuler->func) {VecDestroy(beuler->func);}
202: if (beuler->rhs) {VecDestroy(beuler->rhs);}
203: PetscFree(beuler);
204: return(0);
205: }
207: /*
208: This defines the nonlinear equation that is to be solved with SNES
209: 1/dt* (U^{n+1} - U^{n}) - F(U^{n+1})
210: */
213: PetscErrorCode TSBEulerFunction(SNES snes,Vec x,Vec y,void *ctx)
214: {
215: TS ts = (TS) ctx;
216: PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1;
218: PetscInt i,n;
221: /* apply user-provided function */
222: TSComputeRHSFunction(ts,ts->ptime,x,y);
223: VecGetArray(ts->vec_sol,&un);
224: VecGetArray(x,&unp1);
225: VecGetArray(y,&Funp1);
226: VecGetLocalSize(x,&n);
227: for (i=0; i<n; i++) {
228: Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i];
229: }
230: VecRestoreArray(ts->vec_sol,&un);
231: VecRestoreArray(x,&unp1);
232: VecRestoreArray(y,&Funp1);
233: return(0);
234: }
236: /*
237: This constructs the Jacobian needed for SNES
238: J = I/dt - J_{F} where J_{F} is the given Jacobian of F at t_{n+1}.
239: x - input vector
240: AA - Jacobian matrix
241: BB - preconditioner matrix, usually the same as AA
242: */
245: PetscErrorCode TSBEulerJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
246: {
247: TS ts = (TS) ctx;
251: /* construct user's Jacobian */
252: TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);
254: /* shift and scale Jacobian */
255: /* this test is a undesirable hack, we assume that if it is MATMFFD then it is
256: obtained from -snes_mf_operator and there is computed directly from the
257: FormFunction() SNES is given and therefor does not need to be shifted/scaled
258: BUT maybe it could be MATMFFD and does require shift in some other case??? */
259: TSSetKSPOperators_BEuler(ts);
260: return(0);
261: }
263: /* ------------------------------------------------------------*/
266: static PetscErrorCode TSSetUp_BEuler_Linear_Constant_Matrix(TS ts)
267: {
268: TS_BEuler *beuler = (TS_BEuler*)ts->data;
272: VecDuplicate(ts->vec_sol,&beuler->update);
273: VecDuplicate(ts->vec_sol,&beuler->rhs);
274:
275: /* build linear system to be solved - should move into TSStep() if dt changes! */
276: /* Set ts->A = ts->Arhs = 1/dt*Alhs - Arhs, used in KSPSolve() */
277: TSSetKSPOperators_BEuler(ts);
278: KSPSetOperators(ts->ksp,ts->A,ts->A,SAME_NONZERO_PATTERN);
279: return(0);
280: }
284: static PetscErrorCode TSSetUp_BEuler_Linear_Variable_Matrix(TS ts)
285: {
286: TS_BEuler *beuler = (TS_BEuler*)ts->data;
290: VecDuplicate(ts->vec_sol,&beuler->update);
291: VecDuplicate(ts->vec_sol,&beuler->rhs);
292: return(0);
293: }
297: static PetscErrorCode TSSetUp_BEuler_Nonlinear(TS ts)
298: {
299: TS_BEuler *beuler = (TS_BEuler*)ts->data;
303: VecDuplicate(ts->vec_sol,&beuler->update);
304: VecDuplicate(ts->vec_sol,&beuler->func);
305: SNESSetFunction(ts->snes,beuler->func,TSBEulerFunction,ts);
306: SNESSetJacobian(ts->snes,ts->Arhs,ts->B,TSBEulerJacobian,ts);
307: return(0);
308: }
309: /*------------------------------------------------------------*/
313: static PetscErrorCode TSSetFromOptions_BEuler_Linear(TS ts)
314: {
316: return(0);
317: }
321: static PetscErrorCode TSSetFromOptions_BEuler_Nonlinear(TS ts)
322: {
324: return(0);
325: }
329: static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
330: {
332: return(0);
333: }
335: /* ------------------------------------------------------------ */
336: /*MC
337: TSBEULER - ODE solver using the implicit backward Euler method
339: Level: beginner
341: .seealso: TSCreate(), TS, TSSetType(), TSEULER
343: M*/
347: PetscErrorCode TSCreate_BEuler(TS ts)
348: {
349: TS_BEuler *beuler;
353: ts->ops->destroy = TSDestroy_BEuler;
354: ts->ops->view = TSView_BEuler;
356: if (ts->problem_type == TS_LINEAR) {
357: if (!ts->Arhs) {
358: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set rhs matrix for linear problem");
359: }
360: if (!ts->ops->rhsmatrix) {
361: ts->ops->setup = TSSetUp_BEuler_Linear_Constant_Matrix;
362: ts->ops->step = TSStep_BEuler_Linear_Constant_Matrix;
363: } else {
364: ts->ops->setup = TSSetUp_BEuler_Linear_Variable_Matrix;
365: ts->ops->step = TSStep_BEuler_Linear_Variable_Matrix;
366: }
367: ts->ops->setfromoptions = TSSetFromOptions_BEuler_Linear;
368: KSPCreate(((PetscObject)ts)->comm,&ts->ksp);
369: PetscObjectIncrementTabLevel((PetscObject)ts->ksp,(PetscObject)ts,1);
370: KSPSetInitialGuessNonzero(ts->ksp,PETSC_TRUE);
371: } else if (ts->problem_type == TS_NONLINEAR) {
372: ts->ops->setup = TSSetUp_BEuler_Nonlinear;
373: ts->ops->step = TSStep_BEuler_Nonlinear;
374: ts->ops->setfromoptions = TSSetFromOptions_BEuler_Nonlinear;
375: SNESCreate(((PetscObject)ts)->comm,&ts->snes);
376: PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);
377: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"No such problem");
379: PetscNewLog(ts,TS_BEuler,&beuler);
380: ts->data = (void*)beuler;
382: return(0);
383: }