Actual source code: ex25.c
petsc-3.13.6 2020-09-29
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);
46: int main(int argc,char **argv)
47: {
48: PetscInt cycle;
51: MPI_Init(&argc,&argv);if (ierr) return ierr;
52: for (cycle=0; cycle<4; cycle++) {
53: Brusselator(argc,argv,cycle);
54: if (ierr) return 1;
55: }
56: MPI_Finalize();
57: return ierr;
58: }
60: PetscErrorCode Brusselator(int argc,char **argv,PetscInt cycle)
61: {
62: TS ts; /* nonlinear solver */
63: Vec X; /* solution, residual vectors */
64: Mat J; /* Jacobian matrix */
65: PetscInt steps,mx;
66: PetscErrorCode ierr;
67: DM da;
68: PetscReal ftime,hx,dt,xmax,xmin;
69: struct _User user; /* user-defined work context */
70: TSConvergedReason reason;
72: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
74: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75: Create distributed array (DMDA) to manage parallel grid and vectors
76: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da);
78: DMSetFromOptions(da);
79: DMSetUp(da);
81: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: Extract global vectors from DMDA;
83: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: DMCreateGlobalVector(da,&X);
86: /* Initialize user Section 1.5 Writing Application Codes with PETSc context */
87: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
88: {
89: user.A = 1;
90: user.B = 3;
91: user.alpha = 0.1;
92: user.uleft = 1;
93: user.uright = 1;
94: user.vleft = 3;
95: user.vright = 3;
96: PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL);
97: PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL);
98: PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL);
99: PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL);
100: PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL);
101: PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL);
102: PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL);
103: }
104: PetscOptionsEnd();
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Create timestepping solver context
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109: TSCreate(PETSC_COMM_WORLD,&ts);
110: TSSetDM(ts,da);
111: TSSetType(ts,TSARKIMEX);
112: TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
113: TSSetIFunction(ts,NULL,FormIFunction,&user);
114: DMSetMatType(da,MATAIJ);
115: DMCreateMatrix(da,&J);
116: TSSetIJacobian(ts,J,J,FormIJacobian,&user);
118: ftime = 1.0;
119: TSSetMaxTime(ts,ftime);
120: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
122: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123: Set initial conditions
124: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125: FormInitialSolution(ts,X,&user);
126: TSSetSolution(ts,X);
127: VecGetSize(X,&mx);
128: hx = 1.0/(PetscReal)(mx/2-1);
129: dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */
130: dt *= PetscPowRealInt(0.2,cycle); /* Shrink the time step in convergence study. */
131: TSSetTimeStep(ts,dt);
132: TSSetTolerances(ts,1e-3*PetscPowRealInt(0.5,cycle),NULL,1e-3*PetscPowRealInt(0.5,cycle),NULL);
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Set runtime options
136: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137: TSSetFromOptions(ts);
139: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140: Solve nonlinear system
141: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142: TSSolve(ts,X);
143: TSGetSolveTime(ts,&ftime);
144: TSGetStepNumber(ts,&steps);
145: TSGetConvergedReason(ts,&reason);
146: VecMin(X,NULL,&xmin);
147: VecMax(X,NULL,&xmax);
148: 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);
150: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151: Free work space.
152: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153: MatDestroy(&J);
154: VecDestroy(&X);
155: TSDestroy(&ts);
156: DMDestroy(&da);
157: PetscFinalize();
158: return ierr;
159: }
161: static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
162: {
163: User user = (User)ptr;
164: DM da;
165: DMDALocalInfo info;
166: PetscInt i;
167: Field *x,*xdot,*f;
168: PetscReal hx;
169: Vec Xloc;
173: TSGetDM(ts,&da);
174: DMDAGetLocalInfo(da,&info);
175: hx = 1.0/(PetscReal)(info.mx-1);
177: /*
178: Scatter ghost points to local vector,using the 2-step process
179: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
180: By placing code between these two statements, computations can be
181: done while messages are in transition.
182: */
183: DMGetLocalVector(da,&Xloc);
184: DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
185: DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);
187: /* Get pointers to vector data */
188: DMDAVecGetArrayRead(da,Xloc,&x);
189: DMDAVecGetArrayRead(da,Xdot,&xdot);
190: DMDAVecGetArray(da,F,&f);
192: /* Compute function over the locally owned part of the grid */
193: for (i=info.xs; i<info.xs+info.xm; i++) {
194: if (i == 0) {
195: f[i].u = hx * (x[i].u - user->uleft);
196: f[i].v = hx * (x[i].v - user->vleft);
197: } else if (i == info.mx-1) {
198: f[i].u = hx * (x[i].u - user->uright);
199: f[i].v = hx * (x[i].v - user->vright);
200: } else {
201: f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx;
202: f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx;
203: }
204: }
206: /* Restore vectors */
207: DMDAVecRestoreArrayRead(da,Xloc,&x);
208: DMDAVecRestoreArrayRead(da,Xdot,&xdot);
209: DMDAVecRestoreArray(da,F,&f);
210: DMRestoreLocalVector(da,&Xloc);
211: return(0);
212: }
214: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
215: {
216: User user = (User)ptr;
217: DM da;
218: DMDALocalInfo info;
219: PetscInt i;
220: PetscReal hx;
221: Field *x,*f;
225: TSGetDM(ts,&da);
226: DMDAGetLocalInfo(da,&info);
227: hx = 1.0/(PetscReal)(info.mx-1);
229: /* Get pointers to vector data */
230: DMDAVecGetArrayRead(da,X,&x);
231: DMDAVecGetArray(da,F,&f);
233: /* Compute function over the locally owned part of the grid */
234: for (i=info.xs; i<info.xs+info.xm; i++) {
235: PetscScalar u = x[i].u,v = x[i].v;
236: f[i].u = hx*(user->A + u*u*v - (user->B+1)*u);
237: f[i].v = hx*(user->B*u - u*u*v);
238: }
240: /* Restore vectors */
241: DMDAVecRestoreArrayRead(da,X,&x);
242: DMDAVecRestoreArray(da,F,&f);
243: return(0);
244: }
246: /* --------------------------------------------------------------------- */
247: /*
248: IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
249: */
250: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
251: {
252: User user = (User)ptr;
254: DMDALocalInfo info;
255: PetscInt i;
256: PetscReal hx;
257: DM da;
258: Field *x,*xdot;
261: TSGetDM(ts,&da);
262: DMDAGetLocalInfo(da,&info);
263: hx = 1.0/(PetscReal)(info.mx-1);
265: /* Get pointers to vector data */
266: DMDAVecGetArrayRead(da,X,&x);
267: DMDAVecGetArrayRead(da,Xdot,&xdot);
269: /* Compute function over the locally owned part of the grid */
270: for (i=info.xs; i<info.xs+info.xm; i++) {
271: if (i == 0 || i == info.mx-1) {
272: const PetscInt row = i,col = i;
273: const PetscScalar vals[2][2] = {{hx,0},{0,hx}};
274: MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES);
275: } else {
276: const PetscInt row = i,col[] = {i-1,i,i+1};
277: const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
278: const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
279: {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
280: MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);
281: }
282: }
284: /* Restore vectors */
285: DMDAVecRestoreArrayRead(da,X,&x);
286: DMDAVecRestoreArrayRead(da,Xdot,&xdot);
288: MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
289: MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
290: if (J != Jpre) {
291: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
292: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
293: }
294: return(0);
295: }
297: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
298: {
299: User user = (User)ctx;
300: DM da;
301: PetscInt i;
302: DMDALocalInfo info;
303: Field *x;
304: PetscReal hx;
308: TSGetDM(ts,&da);
309: DMDAGetLocalInfo(da,&info);
310: hx = 1.0/(PetscReal)(info.mx-1);
312: /* Get pointers to vector data */
313: DMDAVecGetArray(da,X,&x);
315: /* Compute function over the locally owned part of the grid */
316: for (i=info.xs; i<info.xs+info.xm; i++) {
317: PetscReal xi = i*hx;
318: x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi);
319: x[i].v = user->vleft*(1.-xi) + user->vright*xi;
320: }
321: DMDAVecRestoreArray(da,X,&x);
322: return(0);
323: }
325: /*TEST
327: build:
328: requires: c99
330: test:
331: args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3
333: test:
334: suffix: 2
335: args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3
337: TEST*/