Actual source code: ex25.c
petsc-3.5.4 2015-05-23
1: static const char help[] = "Call PetscInitialize multiple times.\n";
2: /*
3: This example is based on the Brusselator tutorial of the same name, but tests multiple calls to PetscInitialize().
4: This is a bad "convergence study" because it only compares min and max values of the solution rather than comparing
5: norms of the errors. For convergence studies, we recommend invoking PetscInitialize() only once and comparing norms
6: of errors (perhaps estimated using an accurate reference solution).
8: Time-dependent Brusselator reaction-diffusion PDE in 1d. Demonstrates IMEX methods and multiple solves.
10: u_t - alpha u_xx = A + u^2 v - (B+1) u
11: v_t - alpha v_xx = B u - u^2 v
12: 0 < x < 1;
13: A = 1, B = 3, alpha = 1/10
15: Initial conditions:
16: u(x,0) = 1 + sin(2 pi x)
17: v(x,0) = 3
19: Boundary conditions:
20: u(0,t) = u(1,t) = 1
21: v(0,t) = v(1,t) = 3
22: */
24: #include <petscdm.h>
25: #include <petscdmda.h>
26: #include <petscts.h>
28: typedef struct {
29: PetscScalar u,v;
30: } Field;
32: typedef struct _User *User;
33: struct _User {
34: PetscReal A,B; /* Reaction coefficients */
35: PetscReal alpha; /* Diffusion coefficient */
36: PetscReal uleft,uright; /* Dirichlet boundary conditions */
37: PetscReal vleft,vright; /* Dirichlet boundary conditions */
38: };
40: static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
41: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
42: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
43: static PetscErrorCode FormInitialSolution(TS,Vec,void*);
44: static int Brusselator(int,char**,PetscInt);
48: int main(int argc,char **argv)
49: {
50: PetscInt cycle;
51: MPI_Init(&argc,&argv);
52: for (cycle=0; cycle<4; cycle++) {
53: Brusselator(argc,argv,cycle);
54: }
55: MPI_Finalize();
56: return 0;
57: }
61: int Brusselator(int argc,char **argv,PetscInt cycle)
62: {
63: TS ts; /* nonlinear solver */
64: Vec X; /* solution, residual vectors */
65: Mat J; /* Jacobian matrix */
66: PetscInt steps,maxsteps,mx;
67: PetscErrorCode ierr;
68: DM da;
69: PetscReal ftime,hx,dt,xmax,xmin;
70: struct _User user; /* user-defined work context */
71: TSConvergedReason reason;
73: PetscInitialize(&argc,&argv,(char*)0,help);
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Create distributed array (DMDA) to manage parallel grid and vectors
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,-11,2,2,NULL,&da);
80: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: Extract global vectors from DMDA;
82: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: DMCreateGlobalVector(da,&X);
85: /* Initialize user application context */
86: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
87: {
88: user.A = 1;
89: user.B = 3;
90: user.alpha = 0.1;
91: user.uleft = 1;
92: user.uright = 1;
93: user.vleft = 3;
94: user.vright = 3;
95: PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL);
96: PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL);
97: PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL);
98: PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL);
99: PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL);
100: PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL);
101: PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL);
102: }
103: PetscOptionsEnd();
105: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: Create timestepping solver context
107: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: TSCreate(PETSC_COMM_WORLD,&ts);
109: TSSetDM(ts,da);
110: TSSetType(ts,TSARKIMEX);
111: TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
112: TSSetIFunction(ts,NULL,FormIFunction,&user);
113: DMSetMatType(da,MATAIJ);
114: DMCreateMatrix(da,&J);
115: TSSetIJacobian(ts,J,J,FormIJacobian,&user);
117: ftime = 1.0;
118: maxsteps = 10000;
119: TSSetDuration(ts,maxsteps,ftime);
121: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122: Set initial conditions
123: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124: FormInitialSolution(ts,X,&user);
125: TSSetSolution(ts,X);
126: VecGetSize(X,&mx);
127: hx = 1.0/(PetscReal)(mx/2-1);
128: dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */
129: dt *= PetscPowRealInt(0.2,cycle); /* Shrink the time step in convergence study. */
130: TSSetInitialTimeStep(ts,0.0,dt);
131: TSSetTolerances(ts,1e-3*PetscPowRealInt(0.5,cycle),NULL,1e-3*PetscPowRealInt(0.5,cycle),NULL);
133: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134: Set runtime options
135: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136: TSSetFromOptions(ts);
138: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139: Solve nonlinear system
140: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141: TSSolve(ts,X);
142: TSGetSolveTime(ts,&ftime);
143: TSGetTimeStepNumber(ts,&steps);
144: TSGetConvergedReason(ts,&reason);
145: VecMin(X,NULL,&xmin);
146: VecMax(X,NULL,&xmax);
147: PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after % 3D steps. Range [%6.4f,%6.4f]\n",TSConvergedReasons[reason],(double)ftime,steps,(double)xmin,(double)xmax);
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Free work space.
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152: MatDestroy(&J);
153: VecDestroy(&X);
154: TSDestroy(&ts);
155: DMDestroy(&da);
156: PetscFinalize();
157: return 0;
158: }
162: static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
163: {
164: User user = (User)ptr;
165: DM da;
166: DMDALocalInfo info;
167: PetscInt i;
168: Field *x,*xdot,*f;
169: PetscReal hx;
170: Vec Xloc;
174: TSGetDM(ts,&da);
175: DMDAGetLocalInfo(da,&info);
176: hx = 1.0/(PetscReal)(info.mx-1);
178: /*
179: Scatter ghost points to local vector,using the 2-step process
180: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
181: By placing code between these two statements, computations can be
182: done while messages are in transition.
183: */
184: DMGetLocalVector(da,&Xloc);
185: DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
186: DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);
188: /* Get pointers to vector data */
189: DMDAVecGetArray(da,Xloc,&x);
190: DMDAVecGetArray(da,Xdot,&xdot);
191: DMDAVecGetArray(da,F,&f);
193: /* Compute function over the locally owned part of the grid */
194: for (i=info.xs; i<info.xs+info.xm; i++) {
195: if (i == 0) {
196: f[i].u = hx * (x[i].u - user->uleft);
197: f[i].v = hx * (x[i].v - user->vleft);
198: } else if (i == info.mx-1) {
199: f[i].u = hx * (x[i].u - user->uright);
200: f[i].v = hx * (x[i].v - user->vright);
201: } else {
202: f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx;
203: f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx;
204: }
205: }
207: /* Restore vectors */
208: DMDAVecRestoreArray(da,Xloc,&x);
209: DMDAVecRestoreArray(da,Xdot,&xdot);
210: DMDAVecRestoreArray(da,F,&f);
211: DMRestoreLocalVector(da,&Xloc);
212: return(0);
213: }
217: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
218: {
219: User user = (User)ptr;
220: DM da;
221: DMDALocalInfo info;
222: PetscInt i;
223: PetscReal hx;
224: Field *x,*f;
228: TSGetDM(ts,&da);
229: DMDAGetLocalInfo(da,&info);
230: hx = 1.0/(PetscReal)(info.mx-1);
232: /* Get pointers to vector data */
233: DMDAVecGetArray(da,X,&x);
234: DMDAVecGetArray(da,F,&f);
236: /* Compute function over the locally owned part of the grid */
237: for (i=info.xs; i<info.xs+info.xm; i++) {
238: PetscScalar u = x[i].u,v = x[i].v;
239: f[i].u = hx*(user->A + u*u*v - (user->B+1)*u);
240: f[i].v = hx*(user->B*u - u*u*v);
241: }
243: /* Restore vectors */
244: DMDAVecRestoreArray(da,X,&x);
245: DMDAVecRestoreArray(da,F,&f);
246: return(0);
247: }
249: /* --------------------------------------------------------------------- */
250: /*
251: IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
252: */
255: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
256: {
257: User user = (User)ptr;
259: DMDALocalInfo info;
260: PetscInt i;
261: PetscReal hx;
262: DM da;
263: Field *x,*xdot;
266: TSGetDM(ts,&da);
267: DMDAGetLocalInfo(da,&info);
268: hx = 1.0/(PetscReal)(info.mx-1);
270: /* Get pointers to vector data */
271: DMDAVecGetArray(da,X,&x);
272: DMDAVecGetArray(da,Xdot,&xdot);
274: /* Compute function over the locally owned part of the grid */
275: for (i=info.xs; i<info.xs+info.xm; i++) {
276: if (i == 0 || i == info.mx-1) {
277: const PetscInt row = i,col = i;
278: const PetscScalar vals[2][2] = {{hx,0},{0,hx}};
279: MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES);
280: } else {
281: const PetscInt row = i,col[] = {i-1,i,i+1};
282: const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
283: const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
284: {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
285: MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);
286: }
287: }
289: /* Restore vectors */
290: DMDAVecRestoreArray(da,X,&x);
291: DMDAVecRestoreArray(da,Xdot,&xdot);
293: MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
294: MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
295: if (J != Jpre) {
296: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
297: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
298: }
299: return(0);
300: }
304: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
305: {
306: User user = (User)ctx;
307: DM da;
308: PetscInt i;
309: DMDALocalInfo info;
310: Field *x;
311: PetscReal hx;
315: TSGetDM(ts,&da);
316: DMDAGetLocalInfo(da,&info);
317: hx = 1.0/(PetscReal)(info.mx-1);
319: /* Get pointers to vector data */
320: DMDAVecGetArray(da,X,&x);
322: /* Compute function over the locally owned part of the grid */
323: for (i=info.xs; i<info.xs+info.xm; i++) {
324: PetscReal xi = i*hx;
325: x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi);
326: x[i].v = user->vleft*(1.-xi) + user->vright*xi;
327: }
328: DMDAVecRestoreArray(da,X,&x);
329: return(0);
330: }