Actual source code: ex25.c
petsc-3.6.1 2015-08-06
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 <petscdm.h>
18: #include <petscdmda.h>
19: #include <petscts.h>
21: typedef struct {
22: PetscScalar u,v;
23: } Field;
25: typedef struct _User *User;
26: struct _User {
27: PetscReal A,B; /* Reaction coefficients */
28: PetscReal alpha; /* Diffusion coefficient */
29: PetscReal uleft,uright; /* Dirichlet boundary conditions */
30: PetscReal vleft,vright; /* Dirichlet boundary conditions */
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,void*);
36: static PetscErrorCode FormInitialSolution(TS,Vec,void*);
40: int main(int argc,char **argv)
41: {
42: TS ts; /* nonlinear solver */
43: Vec X; /* solution, residual vectors */
44: Mat J; /* Jacobian matrix */
45: PetscInt steps,maxsteps,mx;
46: PetscErrorCode ierr;
47: DM da;
48: PetscReal ftime,hx,dt;
49: struct _User user; /* user-defined work context */
50: TSConvergedReason reason;
52: PetscInitialize(&argc,&argv,(char*)0,help);
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Create distributed array (DMDA) to manage parallel grid and vectors
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,-11,2,2,NULL,&da);
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Extract global vectors from DMDA;
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62: DMCreateGlobalVector(da,&X);
64: /* Initialize user application context */
65: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
66: {
67: user.A = 1;
68: user.B = 3;
69: user.alpha = 0.02;
70: user.uleft = 1;
71: user.uright = 1;
72: user.vleft = 3;
73: user.vright = 3;
74: PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL);
75: PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL);
76: PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL);
77: PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL);
78: PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL);
79: PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL);
80: PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL);
81: }
82: PetscOptionsEnd();
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Create timestepping solver context
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: TSCreate(PETSC_COMM_WORLD,&ts);
88: TSSetDM(ts,da);
89: TSSetType(ts,TSARKIMEX);
90: TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);
91: TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
92: TSSetIFunction(ts,NULL,FormIFunction,&user);
93: DMSetMatType(da,MATAIJ);
94: DMCreateMatrix(da,&J);
95: TSSetIJacobian(ts,J,J,FormIJacobian,&user);
97: ftime = 10.0;
98: maxsteps = 10000;
99: TSSetDuration(ts,maxsteps,ftime);
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Set initial conditions
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: FormInitialSolution(ts,X,&user);
105: TSSetSolution(ts,X);
106: VecGetSize(X,&mx);
107: hx = 1.0/(PetscReal)(mx/2-1);
108: dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */
109: TSSetInitialTimeStep(ts,0.0,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: TSGetTimeStepNumber(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 0;
134: }
138: static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
139: {
140: User user = (User)ptr;
141: DM da;
142: DMDALocalInfo info;
143: PetscInt i;
144: Field *x,*xdot,*f;
145: PetscReal hx;
146: Vec Xloc;
150: TSGetDM(ts,&da);
151: DMDAGetLocalInfo(da,&info);
152: hx = 1.0/(PetscReal)(info.mx-1);
154: /*
155: Scatter ghost points to local vector,using the 2-step process
156: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
157: By placing code between these two statements, computations can be
158: done while messages are in transition.
159: */
160: DMGetLocalVector(da,&Xloc);
161: DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
162: DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);
164: /* Get pointers to vector data */
165: DMDAVecGetArrayRead(da,Xloc,&x);
166: DMDAVecGetArrayRead(da,Xdot,&xdot);
167: DMDAVecGetArray(da,F,&f);
169: /* Compute function over the locally owned part of the grid */
170: for (i=info.xs; i<info.xs+info.xm; i++) {
171: if (i == 0) {
172: f[i].u = hx * (x[i].u - user->uleft);
173: f[i].v = hx * (x[i].v - user->vleft);
174: } else if (i == info.mx-1) {
175: f[i].u = hx * (x[i].u - user->uright);
176: f[i].v = hx * (x[i].v - user->vright);
177: } else {
178: f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx;
179: f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx;
180: }
181: }
183: /* Restore vectors */
184: DMDAVecRestoreArrayRead(da,Xloc,&x);
185: DMDAVecRestoreArrayRead(da,Xdot,&xdot);
186: DMDAVecRestoreArray(da,F,&f);
187: DMRestoreLocalVector(da,&Xloc);
188: return(0);
189: }
193: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
194: {
195: User user = (User)ptr;
196: DM da;
197: DMDALocalInfo info;
198: PetscInt i;
199: PetscReal hx;
200: Field *x,*f;
204: TSGetDM(ts,&da);
205: DMDAGetLocalInfo(da,&info);
206: hx = 1.0/(PetscReal)(info.mx-1);
208: /* Get pointers to vector data */
209: DMDAVecGetArrayRead(da,X,&x);
210: DMDAVecGetArray(da,F,&f);
212: /* Compute function over the locally owned part of the grid */
213: for (i=info.xs; i<info.xs+info.xm; i++) {
214: PetscScalar u = x[i].u,v = x[i].v;
215: f[i].u = hx*(user->A + u*u*v - (user->B+1)*u);
216: f[i].v = hx*(user->B*u - u*u*v);
217: }
219: /* Restore vectors */
220: DMDAVecRestoreArrayRead(da,X,&x);
221: DMDAVecRestoreArray(da,F,&f);
222: return(0);
223: }
225: /* --------------------------------------------------------------------- */
226: /*
227: IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
228: */
231: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
232: {
233: User user = (User)ptr;
235: DMDALocalInfo info;
236: PetscInt i;
237: PetscReal hx;
238: DM da;
239: Field *x,*xdot;
242: TSGetDM(ts,&da);
243: DMDAGetLocalInfo(da,&info);
244: hx = 1.0/(PetscReal)(info.mx-1);
246: /* Get pointers to vector data */
247: DMDAVecGetArrayRead(da,X,&x);
248: DMDAVecGetArrayRead(da,Xdot,&xdot);
250: /* Compute function over the locally owned part of the grid */
251: for (i=info.xs; i<info.xs+info.xm; i++) {
252: if (i == 0 || i == info.mx-1) {
253: const PetscInt row = i,col = i;
254: const PetscScalar vals[2][2] = {{hx,0},{0,hx}};
255: MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES);
256: } else {
257: const PetscInt row = i,col[] = {i-1,i,i+1};
258: const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
259: const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
260: {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
261: MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);
262: }
263: }
265: /* Restore vectors */
266: DMDAVecRestoreArrayRead(da,X,&x);
267: DMDAVecRestoreArrayRead(da,Xdot,&xdot);
269: MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
270: MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
271: if (J != Jpre) {
272: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
273: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
274: }
275: return(0);
276: }
280: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
281: {
282: User user = (User)ctx;
283: DM da;
284: PetscInt i;
285: DMDALocalInfo info;
286: Field *x;
287: PetscReal hx;
291: TSGetDM(ts,&da);
292: DMDAGetLocalInfo(da,&info);
293: hx = 1.0/(PetscReal)(info.mx-1);
295: /* Get pointers to vector data */
296: DMDAVecGetArray(da,X,&x);
298: /* Compute function over the locally owned part of the grid */
299: for (i=info.xs; i<info.xs+info.xm; i++) {
300: PetscReal xi = i*hx;
301: x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi);
302: x[i].v = user->vleft*(1.-xi) + user->vright*xi;
303: }
304: DMDAVecRestoreArray(da,X,&x);
305: return(0);
306: }