Actual source code: ex25.c
petsc-3.4.5 2014-06-29
1: static const char help[] = "Time-dependent Brusselator reaction-diffusion PDE in 1d. Demonstrates IMEX methods.\n";
2: /*
3: u_t - alpha u_xx = A + u^2 v - (B+1) u
4: v_t - alpha v_xx = B u - u^2 v
5: 0 < x < 1;
6: A = 1, B = 3, alpha = 1/50
8: Initial conditions:
9: u(x,0) = 1 + sin(2 pi x)
10: v(x,0) = 3
12: Boundary conditions:
13: u(0,t) = u(1,t) = 1
14: v(0,t) = v(1,t) = 3
15: */
17: #include <petscdmda.h>
18: #include <petscts.h>
20: typedef struct {
21: PetscScalar u,v;
22: } Field;
24: typedef struct _User *User;
25: struct _User {
26: PetscReal A,B; /* Reaction coefficients */
27: PetscReal alpha; /* Diffusion coefficient */
28: PetscReal uleft,uright; /* Dirichlet boundary conditions */
29: PetscReal vleft,vright; /* Dirichlet boundary conditions */
30: };
32: static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
33: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
34: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
35: static PetscErrorCode FormInitialSolution(TS,Vec,void*);
39: int main(int argc,char **argv)
40: {
41: TS ts; /* nonlinear solver */
42: Vec X; /* solution, residual vectors */
43: Mat J; /* Jacobian matrix */
44: PetscInt steps,maxsteps,mx;
45: PetscErrorCode ierr;
46: DM da;
47: PetscReal ftime,hx,dt;
48: struct _User user; /* user-defined work context */
49: TSConvergedReason reason;
51: PetscInitialize(&argc,&argv,(char*)0,help);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Create distributed array (DMDA) to manage parallel grid and vectors
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-11,2,2,NULL,&da);
58: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: Extract global vectors from DMDA;
60: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: DMCreateGlobalVector(da,&X);
63: /* Initialize user application context */
64: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
65: {
66: user.A = 1;
67: user.B = 3;
68: user.alpha = 0.02;
69: user.uleft = 1;
70: user.uright = 1;
71: user.vleft = 3;
72: user.vright = 3;
73: PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL);
74: PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL);
75: PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL);
76: PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL);
77: PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL);
78: PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL);
79: PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL);
80: }
81: PetscOptionsEnd();
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Create timestepping solver context
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: TSCreate(PETSC_COMM_WORLD,&ts);
87: TSSetDM(ts,da);
88: TSSetType(ts,TSARKIMEX);
89: TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
90: TSSetIFunction(ts,NULL,FormIFunction,&user);
91: DMCreateMatrix(da,MATAIJ,&J);
92: TSSetIJacobian(ts,J,J,FormIJacobian,&user);
94: ftime = 10.0;
95: maxsteps = 10000;
96: TSSetDuration(ts,maxsteps,ftime);
98: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99: Set initial conditions
100: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101: FormInitialSolution(ts,X,&user);
102: TSSetSolution(ts,X);
103: VecGetSize(X,&mx);
104: hx = 1.0/(PetscReal)(mx/2-1);
105: dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */
106: TSSetInitialTimeStep(ts,0.0,dt);
108: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109: Set runtime options
110: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111: TSSetFromOptions(ts);
113: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114: Solve nonlinear system
115: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: TSSolve(ts,X);
117: TSGetSolveTime(ts,&ftime);
118: TSGetTimeStepNumber(ts,&steps);
119: TSGetConvergedReason(ts,&reason);
120: PetscPrintf(PETSC_COMM_WORLD,"%s at time %G after %D steps\n",TSConvergedReasons[reason],ftime,steps);
122: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123: Free work space.
124: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125: MatDestroy(&J);
126: VecDestroy(&X);
127: TSDestroy(&ts);
128: DMDestroy(&da);
129: PetscFinalize();
130: return 0;
131: }
135: static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
136: {
137: User user = (User)ptr;
138: DM da;
139: DMDALocalInfo info;
140: PetscInt i;
141: Field *x,*xdot,*f;
142: PetscReal hx;
143: Vec Xloc;
147: TSGetDM(ts,&da);
148: DMDAGetLocalInfo(da,&info);
149: hx = 1.0/(PetscReal)(info.mx-1);
151: /*
152: Scatter ghost points to local vector,using the 2-step process
153: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
154: By placing code between these two statements, computations can be
155: done while messages are in transition.
156: */
157: DMGetLocalVector(da,&Xloc);
158: DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
159: DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);
161: /* Get pointers to vector data */
162: DMDAVecGetArray(da,Xloc,&x);
163: DMDAVecGetArray(da,Xdot,&xdot);
164: DMDAVecGetArray(da,F,&f);
166: /* Compute function over the locally owned part of the grid */
167: for (i=info.xs; i<info.xs+info.xm; i++) {
168: if (i == 0) {
169: f[i].u = hx * (x[i].u - user->uleft);
170: f[i].v = hx * (x[i].v - user->vleft);
171: } else if (i == info.mx-1) {
172: f[i].u = hx * (x[i].u - user->uright);
173: f[i].v = hx * (x[i].v - user->vright);
174: } else {
175: f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx;
176: f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx;
177: }
178: }
180: /* Restore vectors */
181: DMDAVecRestoreArray(da,Xloc,&x);
182: DMDAVecRestoreArray(da,Xdot,&xdot);
183: DMDAVecRestoreArray(da,F,&f);
184: DMRestoreLocalVector(da,&Xloc);
185: return(0);
186: }
190: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
191: {
192: User user = (User)ptr;
193: DM da;
194: DMDALocalInfo info;
195: PetscInt i;
196: PetscReal hx;
197: Field *x,*f;
201: TSGetDM(ts,&da);
202: DMDAGetLocalInfo(da,&info);
203: hx = 1.0/(PetscReal)(info.mx-1);
205: /* Get pointers to vector data */
206: DMDAVecGetArray(da,X,&x);
207: DMDAVecGetArray(da,F,&f);
209: /* Compute function over the locally owned part of the grid */
210: for (i=info.xs; i<info.xs+info.xm; i++) {
211: PetscScalar u = x[i].u,v = x[i].v;
212: f[i].u = hx*(user->A + u*u*v - (user->B+1)*u);
213: f[i].v = hx*(user->B*u - u*u*v);
214: }
216: /* Restore vectors */
217: DMDAVecRestoreArray(da,X,&x);
218: DMDAVecRestoreArray(da,F,&f);
219: return(0);
220: }
222: /* --------------------------------------------------------------------- */
223: /*
224: IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
225: */
228: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat *J,Mat *Jpre,MatStructure *str,void *ptr)
229: {
230: User user = (User)ptr;
232: DMDALocalInfo info;
233: PetscInt i;
234: PetscReal hx;
235: DM da;
236: Field *x,*xdot;
239: TSGetDM(ts,&da);
240: DMDAGetLocalInfo(da,&info);
241: hx = 1.0/(PetscReal)(info.mx-1);
243: /* Get pointers to vector data */
244: DMDAVecGetArray(da,X,&x);
245: DMDAVecGetArray(da,Xdot,&xdot);
247: /* Compute function over the locally owned part of the grid */
248: for (i=info.xs; i<info.xs+info.xm; i++) {
249: if (i == 0 || i == info.mx-1) {
250: const PetscInt row = i,col = i;
251: const PetscScalar vals[2][2] = {{hx,0},{0,hx}};
252: MatSetValuesBlocked(*Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES);
253: } else {
254: const PetscInt row = i,col[] = {i-1,i,i+1};
255: const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
256: const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
257: {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
258: MatSetValuesBlocked(*Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);
259: }
260: }
262: /* Restore vectors */
263: DMDAVecRestoreArray(da,X,&x);
264: DMDAVecRestoreArray(da,Xdot,&xdot);
266: MatAssemblyBegin(*Jpre,MAT_FINAL_ASSEMBLY);
267: MatAssemblyEnd(*Jpre,MAT_FINAL_ASSEMBLY);
268: if (*J != *Jpre) {
269: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
270: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
271: }
272: return(0);
273: }
277: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
278: {
279: User user = (User)ctx;
280: DM da;
281: PetscInt i;
282: DMDALocalInfo info;
283: Field *x;
284: PetscReal hx;
288: TSGetDM(ts,&da);
289: DMDAGetLocalInfo(da,&info);
290: hx = 1.0/(PetscReal)(info.mx-1);
292: /* Get pointers to vector data */
293: DMDAVecGetArray(da,X,&x);
295: /* Compute function over the locally owned part of the grid */
296: for (i=info.xs; i<info.xs+info.xm; i++) {
297: PetscReal xi = i*hx;
298: x[i].u = user->uleft*(1.-xi) + user->uright*xi + sin(2.*PETSC_PI*xi);
299: x[i].v = user->vleft*(1.-xi) + user->vright*xi;
300: }
301: DMDAVecRestoreArray(da,X,&x);
302: return(0);
303: }