Actual source code: ex22.c
petsc-3.13.6 2020-09-29
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
9: Initial conditions:
10: u(x,0) = 1 + s2*x
11: v(x,0) = k0/k1*u(x,0) + s1/k1
13: Upstream boundary conditions:
14: u(0,t) = 1-sin(12*t)^4
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: */
21: #include <petscdm.h>
22: #include <petscdmda.h>
23: #include <petscts.h>
25: typedef PetscScalar Field[2];
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: };
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*);
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;
53: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
55: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56: Create distributed array (DMDA) to manage parallel grid and vectors
57: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da);
59: DMSetFromOptions(da);
60: DMSetUp(da);
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Extract global vectors from DMDA;
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: DMCreateGlobalVector(da,&X);
67: /* Initialize user Section 1.5 Writing Application Codes with PETSc context */
68: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
69: {
70: user.a[0] = 1; PetscOptionsReal("-a0","Advection rate 0","",user.a[0],&user.a[0],NULL);
71: user.a[1] = 0; PetscOptionsReal("-a1","Advection rate 1","",user.a[1],&user.a[1],NULL);
72: user.k[0] = 1e6; PetscOptionsReal("-k0","Reaction rate 0","",user.k[0],&user.k[0],NULL);
73: user.k[1] = 2*user.k[0]; PetscOptionsReal("-k1","Reaction rate 1","",user.k[1],&user.k[1],NULL);
74: user.s[0] = 0; PetscOptionsReal("-s0","Source 0","",user.s[0],&user.s[0],NULL);
75: user.s[1] = 1; PetscOptionsReal("-s1","Source 1","",user.s[1],&user.s[1],NULL);
76: }
77: PetscOptionsEnd();
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Create timestepping solver context
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: TSCreate(PETSC_COMM_WORLD,&ts);
83: TSSetDM(ts,da);
84: TSSetType(ts,TSARKIMEX);
85: TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
86: TSSetIFunction(ts,NULL,FormIFunction,&user);
87: DMSetMatType(da,MATAIJ);
88: DMCreateMatrix(da,&J);
89: TSSetIJacobian(ts,J,J,FormIJacobian,&user);
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: TSGetSNES(ts,&snes);
95: SNESGetLineSearch(snes,&linesearch);
96: SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
98: ftime = .1;
99: TSSetMaxTime(ts,ftime);
100: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Set initial conditions
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: FormInitialSolution(ts,X,&user);
106: TSSetSolution(ts,X);
107: VecGetSize(X,&mx);
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: TSSetTimeStep(ts,dt);
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Set runtime options
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: TSSetFromOptions(ts);
116: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: Solve nonlinear system
118: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119: TSSolve(ts,X);
120: TSGetSolveTime(ts,&ftime);
121: TSGetStepNumber(ts,&steps);
122: TSGetConvergedReason(ts,&reason);
123: PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);
125: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: Free work space.
127: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128: MatDestroy(&J);
129: VecDestroy(&X);
130: TSDestroy(&ts);
131: DMDestroy(&da);
132: PetscFinalize();
133: return ierr;
134: }
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;
147: TSGetDM(ts,&da);
148: DMDAGetLocalInfo(da,&info);
150: /* Get pointers to vector data */
151: DMDAVecGetArrayRead(da,X,(void*)&x);
152: DMDAVecGetArrayRead(da,Xdot,(void*)&xdot);
153: DMDAVecGetArray(da,F,&f);
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: }
161: /* Restore vectors */
162: DMDAVecRestoreArrayRead(da,X,(void*)&x);
163: DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot);
164: DMDAVecRestoreArray(da,F,&f);
165: return(0);
166: }
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;
181: TSGetDM(ts,&da);
182: DMDAGetLocalInfo(da,&info);
183: hx = 1.0/(PetscReal)info.mx;
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: DMGetLocalVector(da,&Xloc);
192: DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
193: DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);
195: /* Get pointers to vector data */
196: DMDAVecGetArrayRead(da,Xloc,(void*)&x);
197: DMDAVecGetArray(da,F,&f);
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: }
214: /* Restore vectors */
215: DMDAVecRestoreArrayRead(da,Xloc,(void*)&x);
216: DMDAVecRestoreArray(da,F,&f);
217: DMRestoreLocalVector(da,&Xloc);
218: return(0);
219: }
221: /* --------------------------------------------------------------------- */
222: /*
223: IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
224: */
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;
229: DMDALocalInfo info;
230: PetscInt i;
231: DM da;
232: const Field *x,*xdot;
235: TSGetDM(ts,&da);
236: DMDAGetLocalInfo(da,&info);
238: /* Get pointers to vector data */
239: DMDAVecGetArrayRead(da,X,(void*)&x);
240: DMDAVecGetArrayRead(da,Xdot,(void*)&xdot);
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];
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: MatSetValuesBlocked(Jpre,1,&i,1,&i,&v[0][0],INSERT_VALUES);
250: }
252: /* Restore vectors */
253: DMDAVecRestoreArrayRead(da,X,(void*)&x);
254: DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot);
256: MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
257: MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
258: if (J != Jpre) {
259: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
260: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
261: }
262: return(0);
263: }
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;
276: TSGetDM(ts,&da);
277: DMDAGetLocalInfo(da,&info);
278: hx = 1.0/(PetscReal)info.mx;
280: /* Get pointers to vector data */
281: DMDAVecGetArray(da,X,&x);
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: DMDAVecRestoreArray(da,X,&x);
290: return(0);
291: }
293: /*TEST
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
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
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
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
315: TEST*/