1 static const char help[] = "Time-dependent advection-reaction PDE in 1d, demonstrates IMEX methods.\n";
2 /*
3 u_t + a1*u_x = -k1*u + k2*v + s1
4 v_t + a2*v_x = k1*u - k2*v + s2
5 0 < x < 1;
6 a1 = 1, k1 = 10^6, s1 = 0,
7 a2 = 0, k2 = 2*k1, s2 = 1
8
9 Initial conditions:
10 u(x,0) = 1 + s2*x
11 v(x,0) = k0/k1*u(x,0) + s1/k1
12
13 Upstream boundary conditions:
14 u(0,t) = 1-sin(12*t)^4
15
16 Useful command-line parameters:
17 -ts_arkimex_fully_implicit - treats advection implicitly, only relevant with -ts_type arkimex (default)
18 -ts_type rosw - use Rosenbrock-W method (linearly implicit IMEX, amortizes assembly and preconditioner setup)
19 */
20
21 #include <petscdm.h>
22 #include <petscdmda.h>
23 #include <petscts.h>
24
25 typedef PetscScalar Field[2];
26
27 typedef struct _User *User;
28 struct _User {
29 PetscReal a[2]; /* Advection speeds */
30 PetscReal k[2]; /* Reaction coefficients */
31 PetscReal s[2]; /* Source terms */
32 };
33
34 static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
35 static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
36 static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
37 static PetscErrorCode FormInitialSolution(TS,Vec,void*);
38
main(int argc,char ** argv)39 int main(int argc,char **argv)
40 {
41 TS ts; /* time integrator */
42 SNES snes; /* nonlinear solver */
43 SNESLineSearch linesearch; /* line search */
44 Vec X; /* solution, residual vectors */
45 Mat J; /* Jacobian matrix */
46 PetscInt steps,mx;
47 PetscErrorCode ierr;
48 DM da;
49 PetscReal ftime,dt;
50 struct _User user; /* user-defined work context */
51 TSConvergedReason reason;
52
53 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
54
55 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56 Create distributed array (DMDA) to manage parallel grid and vectors
57 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58 ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da);CHKERRQ(ierr);
59 ierr = DMSetFromOptions(da);CHKERRQ(ierr);
60 ierr = DMSetUp(da);CHKERRQ(ierr);
61
62 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63 Extract global vectors from DMDA;
64 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65 ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr);
66
67 /* Initialize user application context */
68 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
69 {
70 user.a[0] = 1; ierr = PetscOptionsReal("-a0","Advection rate 0","",user.a[0],&user.a[0],NULL);CHKERRQ(ierr);
71 user.a[1] = 0; ierr = PetscOptionsReal("-a1","Advection rate 1","",user.a[1],&user.a[1],NULL);CHKERRQ(ierr);
72 user.k[0] = 1e6; ierr = PetscOptionsReal("-k0","Reaction rate 0","",user.k[0],&user.k[0],NULL);CHKERRQ(ierr);
73 user.k[1] = 2*user.k[0]; ierr = PetscOptionsReal("-k1","Reaction rate 1","",user.k[1],&user.k[1],NULL);CHKERRQ(ierr);
74 user.s[0] = 0; ierr = PetscOptionsReal("-s0","Source 0","",user.s[0],&user.s[0],NULL);CHKERRQ(ierr);
75 user.s[1] = 1; ierr = PetscOptionsReal("-s1","Source 1","",user.s[1],&user.s[1],NULL);CHKERRQ(ierr);
76 }
77 ierr = PetscOptionsEnd();CHKERRQ(ierr);
78
79 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80 Create timestepping solver context
81 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82 ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
83 ierr = TSSetDM(ts,da);CHKERRQ(ierr);
84 ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr);
85 ierr = TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);CHKERRQ(ierr);
86 ierr = TSSetIFunction(ts,NULL,FormIFunction,&user);CHKERRQ(ierr);
87 ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr);
88 ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr);
89 ierr = TSSetIJacobian(ts,J,J,FormIJacobian,&user);CHKERRQ(ierr);
90
91 /* A line search in the nonlinear solve can fail due to ill-conditioning unless an absolute tolerance is set. Since
92 * this problem is linear, we deactivate the line search. For a linear problem, it is usually recommended to also use
93 * SNESSetType(snes,SNESKSPONLY). */
94 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
95 ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
96 ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
97
98 ftime = .1;
99 ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
100 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
101
102 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103 Set initial conditions
104 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105 ierr = FormInitialSolution(ts,X,&user);CHKERRQ(ierr);
106 ierr = TSSetSolution(ts,X);CHKERRQ(ierr);
107 ierr = VecGetSize(X,&mx);CHKERRQ(ierr);
108 dt = .1 * PetscMax(user.a[0],user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */
109 ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
110
111 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112 Set runtime options
113 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114 ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
115
116 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117 Solve nonlinear system
118 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119 ierr = TSSolve(ts,X);CHKERRQ(ierr);
120 ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
121 ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
122 ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr);
123 ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);CHKERRQ(ierr);
124
125 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126 Free work space.
127 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128 ierr = MatDestroy(&J);CHKERRQ(ierr);
129 ierr = VecDestroy(&X);CHKERRQ(ierr);
130 ierr = TSDestroy(&ts);CHKERRQ(ierr);
131 ierr = DMDestroy(&da);CHKERRQ(ierr);
132 ierr = PetscFinalize();
133 return ierr;
134 }
135
FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void * ptr)136 static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
137 {
138 User user = (User)ptr;
139 DM da;
140 DMDALocalInfo info;
141 PetscInt i;
142 Field *f;
143 const Field *x,*xdot;
144 PetscErrorCode ierr;
145
146 PetscFunctionBeginUser;
147 ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
148 ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
149
150 /* Get pointers to vector data */
151 ierr = DMDAVecGetArrayRead(da,X,(void*)&x);CHKERRQ(ierr);
152 ierr = DMDAVecGetArrayRead(da,Xdot,(void*)&xdot);CHKERRQ(ierr);
153 ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
154
155 /* Compute function over the locally owned part of the grid */
156 for (i=info.xs; i<info.xs+info.xm; i++) {
157 f[i][0] = xdot[i][0] + user->k[0]*x[i][0] - user->k[1]*x[i][1] - user->s[0];
158 f[i][1] = xdot[i][1] - user->k[0]*x[i][0] + user->k[1]*x[i][1] - user->s[1];
159 }
160
161 /* Restore vectors */
162 ierr = DMDAVecRestoreArrayRead(da,X,(void*)&x);CHKERRQ(ierr);
163 ierr = DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot);CHKERRQ(ierr);
164 ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
165 PetscFunctionReturn(0);
166 }
167
FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void * ptr)168 static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
169 {
170 User user = (User)ptr;
171 DM da;
172 Vec Xloc;
173 DMDALocalInfo info;
174 PetscInt i,j;
175 PetscReal hx;
176 Field *f;
177 const Field *x;
178 PetscErrorCode ierr;
179
180 PetscFunctionBeginUser;
181 ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
182 ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
183 hx = 1.0/(PetscReal)info.mx;
184
185 /*
186 Scatter ghost points to local vector,using the 2-step process
187 DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
188 By placing code between these two statements, computations can be
189 done while messages are in transition.
190 */
191 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
192 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
193 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
194
195 /* Get pointers to vector data */
196 ierr = DMDAVecGetArrayRead(da,Xloc,(void*)&x);CHKERRQ(ierr);
197 ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
198
199 /* Compute function over the locally owned part of the grid */
200 for (i=info.xs; i<info.xs+info.xm; i++) {
201 const PetscReal *a = user->a;
202 PetscReal u0t[2];
203 u0t[0] = 1.0 - PetscPowRealInt(PetscSinReal(12*t),4);
204 u0t[1] = 0.0;
205 for (j=0; j<2; j++) {
206 if (i == 0) f[i][j] = a[j]/hx*(1./3*u0t[j] + 0.5*x[i][j] - x[i+1][j] + 1./6*x[i+2][j]);
207 else if (i == 1) f[i][j] = a[j]/hx*(-1./12*u0t[j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]);
208 else if (i == info.mx-2) f[i][j] = a[j]/hx*(-1./6*x[i-2][j] + x[i-1][j] - 0.5*x[i][j] - 1./3*x[i+1][j]);
209 else if (i == info.mx-1) f[i][j] = a[j]/hx*(-x[i][j] + x[i-1][j]);
210 else f[i][j] = a[j]/hx*(-1./12*x[i-2][j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]);
211 }
212 }
213
214 /* Restore vectors */
215 ierr = DMDAVecRestoreArrayRead(da,Xloc,(void*)&x);CHKERRQ(ierr);
216 ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
217 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
218 PetscFunctionReturn(0);
219 }
220
221 /* --------------------------------------------------------------------- */
222 /*
223 IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
224 */
FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void * ptr)225 PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
226 {
227 User user = (User)ptr;
228 PetscErrorCode ierr;
229 DMDALocalInfo info;
230 PetscInt i;
231 DM da;
232 const Field *x,*xdot;
233
234 PetscFunctionBeginUser;
235 ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
236 ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
237
238 /* Get pointers to vector data */
239 ierr = DMDAVecGetArrayRead(da,X,(void*)&x);CHKERRQ(ierr);
240 ierr = DMDAVecGetArrayRead(da,Xdot,(void*)&xdot);CHKERRQ(ierr);
241
242 /* Compute function over the locally owned part of the grid */
243 for (i=info.xs; i<info.xs+info.xm; i++) {
244 const PetscReal *k = user->k;
245 PetscScalar v[2][2];
246
247 v[0][0] = a + k[0]; v[0][1] = -k[1];
248 v[1][0] = -k[0]; v[1][1] = a+k[1];
249 ierr = MatSetValuesBlocked(Jpre,1,&i,1,&i,&v[0][0],INSERT_VALUES);CHKERRQ(ierr);
250 }
251
252 /* Restore vectors */
253 ierr = DMDAVecRestoreArrayRead(da,X,(void*)&x);CHKERRQ(ierr);
254 ierr = DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot);CHKERRQ(ierr);
255
256 ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
257 ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258 if (J != Jpre) {
259 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
260 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
261 }
262 PetscFunctionReturn(0);
263 }
264
FormInitialSolution(TS ts,Vec X,void * ctx)265 PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
266 {
267 User user = (User)ctx;
268 DM da;
269 PetscInt i;
270 DMDALocalInfo info;
271 Field *x;
272 PetscReal hx;
273 PetscErrorCode ierr;
274
275 PetscFunctionBeginUser;
276 ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
277 ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
278 hx = 1.0/(PetscReal)info.mx;
279
280 /* Get pointers to vector data */
281 ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
282
283 /* Compute function over the locally owned part of the grid */
284 for (i=info.xs; i<info.xs+info.xm; i++) {
285 PetscReal r = (i+1)*hx,ik = user->k[1] != 0.0 ? 1.0/user->k[1] : 1.0;
286 x[i][0] = 1 + user->s[1]*r;
287 x[i][1] = user->k[0]*ik*x[i][0] + user->s[1]*ik;
288 }
289 ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
290 PetscFunctionReturn(0);
291 }
292
293 /*TEST
294
295 test:
296 args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_arkimex_type 4 -ts_adapt_type none -ts_dt .005 -ts_max_time .1
297 requires: !single
298
299 test:
300 suffix: 2
301 args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_dt 1e-3 -ts_adapt_type none -ts_dt .005 -ts_max_time .1
302 nsize: 2
303
304 test:
305 suffix: 3
306 args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_rosw_type ra34pw2 -ts_dt 5e-3 -ts_max_time .1 -ts_adapt_type none
307 nsize: 2
308
309 test:
310 suffix: 4
311 args: -ts_type eimex -da_grid_x 200 -ts_eimex_order_adapt true -ts_dt 0.001 -ts_monitor -ts_max_steps 100
312 filter: sed "s/ITS/TIME/g"
313 nsize: 2
314
315 TEST*/
316