Actual source code: ex25.c

petsc-dev 2014-02-02
Report Typos and Errors
  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 <petscdmda.h>
 25: #include <petscts.h>

 27: typedef struct {
 28:   PetscScalar u,v;
 29: } Field;

 31: typedef struct _User *User;
 32: struct _User {
 33:   PetscReal A,B;                /* Reaction coefficients */
 34:   PetscReal alpha;              /* Diffusion coefficient */
 35:   PetscReal uleft,uright;       /* Dirichlet boundary conditions */
 36:   PetscReal vleft,vright;       /* Dirichlet boundary conditions */
 37: };

 39: static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
 40: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
 41: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
 42: static PetscErrorCode FormInitialSolution(TS,Vec,void*);
 43: static int Brusselator(int,char**,PetscInt);

 47: int main(int argc,char **argv)
 48: {
 49:   PetscInt cycle;
 50:   MPI_Init(&argc,&argv);
 51:   for (cycle=0; cycle<4; cycle++) {
 52:     Brusselator(argc,argv,cycle);
 53:   }
 54:   MPI_Finalize();
 55:   return 0;
 56: }

 60: int 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,maxsteps,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);

 74:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75:      Create distributed array (DMDA) to manage parallel grid and vectors
 76:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 77:   DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-11,2,2,NULL,&da);

 79:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:      Extract global vectors from DMDA;
 81:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 82:   DMCreateGlobalVector(da,&X);

 84:   /* Initialize user application context */
 85:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
 86:   {
 87:     user.A      = 1;
 88:     user.B      = 3;
 89:     user.alpha  = 0.1;
 90:     user.uleft  = 1;
 91:     user.uright = 1;
 92:     user.vleft  = 3;
 93:     user.vright = 3;
 94:     PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL);
 95:     PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL);
 96:     PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL);
 97:     PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL);
 98:     PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL);
 99:     PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL);
100:     PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL);
101:   }
102:   PetscOptionsEnd();

104:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105:      Create timestepping solver context
106:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107:   TSCreate(PETSC_COMM_WORLD,&ts);
108:   TSSetDM(ts,da);
109:   TSSetType(ts,TSARKIMEX);
110:   TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
111:   TSSetIFunction(ts,NULL,FormIFunction,&user);
112:   DMSetMatType(da,MATAIJ);
113:   DMCreateMatrix(da,&J);
114:   TSSetIJacobian(ts,J,J,FormIJacobian,&user);

116:   ftime    = 1.0;
117:   maxsteps = 10000;
118:   TSSetDuration(ts,maxsteps,ftime);

120:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121:      Set initial conditions
122:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123:   FormInitialSolution(ts,X,&user);
124:   TSSetSolution(ts,X);
125:   VecGetSize(X,&mx);
126:   hx = 1.0/(PetscReal)(mx/2-1);
127:   dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */
128:   dt *= PetscPowRealInt(0.2,cycle);     /* Shrink the time step in convergence study. */
129:   TSSetInitialTimeStep(ts,0.0,dt);
130:   TSSetTolerances(ts,1e-3*PetscPowRealInt(0.5,cycle),NULL,1e-3*PetscPowRealInt(0.5,cycle),NULL);

132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:      Set runtime options
134:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135:   TSSetFromOptions(ts);

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:      Solve nonlinear system
139:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:   TSSolve(ts,X);
141:   TSGetSolveTime(ts,&ftime);
142:   TSGetTimeStepNumber(ts,&steps);
143:   TSGetConvergedReason(ts,&reason);
144:   VecMin(X,NULL,&xmin);
145:   VecMax(X,NULL,&xmax);
146:   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);

148:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149:      Free work space.
150:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151:   MatDestroy(&J);
152:   VecDestroy(&X);
153:   TSDestroy(&ts);
154:   DMDestroy(&da);
155:   PetscFinalize();
156:   return 0;
157: }

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:   DMDAVecGetArray(da,Xloc,&x);
189:   DMDAVecGetArray(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:   DMDAVecRestoreArray(da,Xloc,&x);
208:   DMDAVecRestoreArray(da,Xdot,&xdot);
209:   DMDAVecRestoreArray(da,F,&f);
210:   DMRestoreLocalVector(da,&Xloc);
211:   return(0);
212: }

216: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
217: {
218:   User           user = (User)ptr;
219:   DM             da;
220:   DMDALocalInfo  info;
221:   PetscInt       i;
222:   PetscReal      hx;
223:   Field          *x,*f;

227:   TSGetDM(ts,&da);
228:   DMDAGetLocalInfo(da,&info);
229:   hx   = 1.0/(PetscReal)(info.mx-1);

231:   /* Get pointers to vector data */
232:   DMDAVecGetArray(da,X,&x);
233:   DMDAVecGetArray(da,F,&f);

235:   /* Compute function over the locally owned part of the grid */
236:   for (i=info.xs; i<info.xs+info.xm; i++) {
237:     PetscScalar u = x[i].u,v = x[i].v;
238:     f[i].u = hx*(user->A + u*u*v - (user->B+1)*u);
239:     f[i].v = hx*(user->B*u - u*u*v);
240:   }

242:   /* Restore vectors */
243:   DMDAVecRestoreArray(da,X,&x);
244:   DMDAVecRestoreArray(da,F,&f);
245:   return(0);
246: }

248: /* --------------------------------------------------------------------- */
249: /*
250:   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
251: */
254: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat *J,Mat *Jpre,MatStructure *str,void *ptr)
255: {
256:   User           user = (User)ptr;
258:   DMDALocalInfo  info;
259:   PetscInt       i;
260:   PetscReal      hx;
261:   DM             da;
262:   Field          *x,*xdot;

265:   TSGetDM(ts,&da);
266:   DMDAGetLocalInfo(da,&info);
267:   hx   = 1.0/(PetscReal)(info.mx-1);

269:   /* Get pointers to vector data */
270:   DMDAVecGetArray(da,X,&x);
271:   DMDAVecGetArray(da,Xdot,&xdot);

273:   /* Compute function over the locally owned part of the grid */
274:   for (i=info.xs; i<info.xs+info.xm; i++) {
275:     if (i == 0 || i == info.mx-1) {
276:       const PetscInt    row        = i,col = i;
277:       const PetscScalar vals[2][2] = {{hx,0},{0,hx}};
278:       MatSetValuesBlocked(*Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES);
279:     } else {
280:       const PetscInt    row           = i,col[] = {i-1,i,i+1};
281:       const PetscScalar dxxL          = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
282:       const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
283:                                          {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
284:       MatSetValuesBlocked(*Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);
285:     }
286:   }

288:   /* Restore vectors */
289:   DMDAVecRestoreArray(da,X,&x);
290:   DMDAVecRestoreArray(da,Xdot,&xdot);

292:   MatAssemblyBegin(*Jpre,MAT_FINAL_ASSEMBLY);
293:   MatAssemblyEnd(*Jpre,MAT_FINAL_ASSEMBLY);
294:   if (*J != *Jpre) {
295:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
296:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
297:   }
298:   return(0);
299: }

303: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
304: {
305:   User           user = (User)ctx;
306:   DM             da;
307:   PetscInt       i;
308:   DMDALocalInfo  info;
309:   Field          *x;
310:   PetscReal      hx;

314:   TSGetDM(ts,&da);
315:   DMDAGetLocalInfo(da,&info);
316:   hx   = 1.0/(PetscReal)(info.mx-1);

318:   /* Get pointers to vector data */
319:   DMDAVecGetArray(da,X,&x);

321:   /* Compute function over the locally owned part of the grid */
322:   for (i=info.xs; i<info.xs+info.xm; i++) {
323:     PetscReal xi = i*hx;
324:     x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi);
325:     x[i].v = user->vleft*(1.-xi) + user->vright*xi;
326:   }
327:   DMDAVecRestoreArray(da,X,&x);
328:   return(0);
329: }