Actual source code: ex22.c
petsc-3.4.5 2014-06-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 <petscdmda.h>
22: #include <petscts.h>
24: typedef PetscScalar Field[2];
26: typedef struct _User *User;
27: struct _User {
28: PetscReal a[2]; /* Advection speeds */
29: PetscReal k[2]; /* Reaction coefficients */
30: PetscReal s[2]; /* Source terms */
31: };
33: static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
34: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
35: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
36: static PetscErrorCode FormInitialSolution(TS,Vec,void*);
40: int main(int argc,char **argv)
41: {
42: TS ts; /* time integrator */
43: SNES snes; /* nonlinear solver */
44: SNESLineSearch linesearch; /* line search */
45: Vec X; /* solution, residual vectors */
46: Mat J; /* Jacobian matrix */
47: PetscInt steps,maxsteps,mx;
48: PetscErrorCode ierr;
49: DM da;
50: PetscReal ftime,dt;
51: struct _User user; /* user-defined work context */
52: TSConvergedReason reason;
54: PetscInitialize(&argc,&argv,(char*)0,help);
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Create distributed array (DMDA) to manage parallel grid and vectors
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-11,2,2,NULL,&da);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Extract global vectors from DMDA;
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: DMCreateGlobalVector(da,&X);
66: /* Initialize user application context */
67: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
68: {
69: user.a[0] = 1; PetscOptionsReal("-a0","Advection rate 0","",user.a[0],&user.a[0],NULL);
70: user.a[1] = 0; PetscOptionsReal("-a1","Advection rate 1","",user.a[1],&user.a[1],NULL);
71: user.k[0] = 1e6; PetscOptionsReal("-k0","Reaction rate 0","",user.k[0],&user.k[0],NULL);
72: user.k[1] = 2*user.k[0]; PetscOptionsReal("-k1","Reaction rate 1","",user.k[1],&user.k[1],NULL);
73: user.s[0] = 0; PetscOptionsReal("-s0","Source 0","",user.s[0],&user.s[0],NULL);
74: user.s[1] = 1; PetscOptionsReal("-s1","Source 1","",user.s[1],&user.s[1],NULL);
75: }
76: PetscOptionsEnd();
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Create timestepping solver context
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: TSCreate(PETSC_COMM_WORLD,&ts);
82: TSSetDM(ts,da);
83: TSSetType(ts,TSARKIMEX);
84: TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
85: TSSetIFunction(ts,NULL,FormIFunction,&user);
86: DMCreateMatrix(da,MATAIJ,&J);
87: TSSetIJacobian(ts,J,J,FormIJacobian,&user);
89: /* A line search in the nonlinear solve can fail due to ill-conditioning unless an absolute tolerance is set. Since
90: * this problem is linear, we deactivate the line search. For a linear problem, it is usually recommended to also use
91: * SNESSetType(snes,SNESKSPONLY). */
92: TSGetSNES(ts,&snes);
93: SNESGetLineSearch(snes,&linesearch);
94: SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
96: ftime = 1.0;
97: maxsteps = 10000;
98: TSSetDuration(ts,maxsteps,ftime);
100: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: Set initial conditions
102: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103: FormInitialSolution(ts,X,&user);
104: TSSetSolution(ts,X);
105: VecGetSize(X,&mx);
106: dt = .1 * PetscMax(user.a[0],user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */
107: TSSetInitialTimeStep(ts,0.0,dt);
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Set runtime options
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112: TSSetFromOptions(ts);
114: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: Solve nonlinear system
116: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117: TSSolve(ts,X);
118: TSGetSolveTime(ts,&ftime);
119: TSGetTimeStepNumber(ts,&steps);
120: TSGetConvergedReason(ts,&reason);
121: PetscPrintf(PETSC_COMM_WORLD,"%s at time %G after %D steps\n",TSConvergedReasons[reason],ftime,steps);
123: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124: Free work space.
125: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126: MatDestroy(&J);
127: VecDestroy(&X);
128: TSDestroy(&ts);
129: DMDestroy(&da);
130: PetscFinalize();
131: return 0;
132: }
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 *x,*xdot,*f;
146: TSGetDM(ts,&da);
147: DMDAGetLocalInfo(da,&info);
149: /* Get pointers to vector data */
150: DMDAVecGetArray(da,X,&x);
151: DMDAVecGetArray(da,Xdot,&xdot);
152: DMDAVecGetArray(da,F,&f);
154: /* Compute function over the locally owned part of the grid */
155: for (i=info.xs; i<info.xs+info.xm; i++) {
156: f[i][0] = xdot[i][0] + user->k[0]*x[i][0] - user->k[1]*x[i][1] - user->s[0];
157: f[i][1] = xdot[i][1] - user->k[0]*x[i][0] + user->k[1]*x[i][1] - user->s[1];
158: }
160: /* Restore vectors */
161: DMDAVecRestoreArray(da,X,&x);
162: DMDAVecRestoreArray(da,Xdot,&xdot);
163: DMDAVecRestoreArray(da,F,&f);
164: return(0);
165: }
169: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
170: {
171: User user = (User)ptr;
172: DM da;
173: Vec Xloc;
174: DMDALocalInfo info;
175: PetscInt i,j;
176: PetscReal hx;
177: Field *x,*f;
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: DMDAVecGetArray(da,Xloc,&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] = {1. - PetscPowScalar(sin(12*t),4.),0};
203: for (j=0; j<2; j++) {
204: 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]);
205: 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]);
206: 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]);
207: else if (i == info.mx-1) f[i][j] = a[j]/hx*(-x[i][j] + x[i-1][j]);
208: 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]);
209: }
210: }
212: /* Restore vectors */
213: DMDAVecRestoreArray(da,Xloc,&x);
214: DMDAVecRestoreArray(da,F,&f);
215: DMRestoreLocalVector(da,&Xloc);
216: return(0);
217: }
219: /* --------------------------------------------------------------------- */
220: /*
221: IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
222: */
225: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat *J,Mat *Jpre,MatStructure *str,void *ptr)
226: {
227: User user = (User)ptr;
229: DMDALocalInfo info;
230: PetscInt i;
231: DM da;
232: Field *x,*xdot;
235: TSGetDM(ts,&da);
236: DMDAGetLocalInfo(da,&info);
238: /* Get pointers to vector data */
239: DMDAVecGetArray(da,X,&x);
240: DMDAVecGetArray(da,Xdot,&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: DMDAVecRestoreArray(da,X,&x);
254: DMDAVecRestoreArray(da,Xdot,&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: }
267: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
268: {
269: User user = (User)ctx;
270: DM da;
271: PetscInt i;
272: DMDALocalInfo info;
273: Field *x;
274: PetscReal hx;
278: TSGetDM(ts,&da);
279: DMDAGetLocalInfo(da,&info);
280: hx = 1.0/(PetscReal)info.mx;
282: /* Get pointers to vector data */
283: DMDAVecGetArray(da,X,&x);
285: /* Compute function over the locally owned part of the grid */
286: for (i=info.xs; i<info.xs+info.xm; i++) {
287: PetscReal r = (i+1)*hx,ik = user->k[1] != 0.0 ? 1.0/user->k[1] : 1.0;
288: x[i][0] = 1 + user->s[1]*r;
289: x[i][1] = user->k[0]*ik*x[i][0] + user->s[1]*ik;
290: }
291: DMDAVecRestoreArray(da,X,&x);
292: return(0);
293: }