Actual source code: ex25.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  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*);

 38: int main(int argc,char **argv)
 39: {
 40:   TS                ts;         /* nonlinear solver */
 41:   Vec               X;          /* solution, residual vectors */
 42:   Mat               J;          /* Jacobian matrix */
 43:   PetscInt          steps,mx;
 44:   PetscErrorCode    ierr;
 45:   DM                da;
 46:   PetscReal         ftime,hx,dt;
 47:   struct _User      user;       /* user-defined work context */
 48:   TSConvergedReason reason;

 50:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 51:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 52:      Create distributed array (DMDA) to manage parallel grid and vectors
 53:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 54:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da);
 55:   DMSetFromOptions(da);
 56:   DMSetUp(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:   TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);
 90:   TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
 91:   TSSetIFunction(ts,NULL,FormIFunction,&user);
 92:   DMSetMatType(da,MATAIJ);
 93:   DMCreateMatrix(da,&J);
 94:   TSSetIJacobian(ts,J,J,FormIJacobian,&user);

 96:   ftime = 10.0;
 97:   TSSetMaxTime(ts,ftime);
 98:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:      Set initial conditions
102:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103:   FormInitialSolution(ts,X,&user);
104:   TSSetSolution(ts,X);
105:   VecGetSize(X,&mx);
106:   hx = 1.0/(PetscReal)(mx/2-1);
107:   dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */
108:   TSSetTimeStep(ts,dt);

110:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111:      Set runtime options
112:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113:   TSSetFromOptions(ts);

115:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116:      Solve nonlinear system
117:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118:   TSSolve(ts,X);
119:   TSGetSolveTime(ts,&ftime);
120:   TSGetStepNumber(ts,&steps);
121:   TSGetConvergedReason(ts,&reason);
122:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);

124:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125:      Free work space.
126:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127:   MatDestroy(&J);
128:   VecDestroy(&X);
129:   TSDestroy(&ts);
130:   DMDestroy(&da);
131:   PetscFinalize();
132:   return ierr;
133: }

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:   DMDAVecGetArrayRead(da,Xloc,&x);
163:   DMDAVecGetArrayRead(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:   DMDAVecRestoreArrayRead(da,Xloc,&x);
182:   DMDAVecRestoreArrayRead(da,Xdot,&xdot);
183:   DMDAVecRestoreArray(da,F,&f);
184:   DMRestoreLocalVector(da,&Xloc);
185:   return(0);
186: }

188: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
189: {
190:   User           user = (User)ptr;
191:   DM             da;
192:   DMDALocalInfo  info;
193:   PetscInt       i;
194:   PetscReal      hx;
195:   Field          *x,*f;

199:   TSGetDM(ts,&da);
200:   DMDAGetLocalInfo(da,&info);
201:   hx   = 1.0/(PetscReal)(info.mx-1);

203:   /* Get pointers to vector data */
204:   DMDAVecGetArrayRead(da,X,&x);
205:   DMDAVecGetArray(da,F,&f);

207:   /* Compute function over the locally owned part of the grid */
208:   for (i=info.xs; i<info.xs+info.xm; i++) {
209:     PetscScalar u = x[i].u,v = x[i].v;
210:     f[i].u = hx*(user->A + u*u*v - (user->B+1)*u);
211:     f[i].v = hx*(user->B*u - u*u*v);
212:   }

214:   /* Restore vectors */
215:   DMDAVecRestoreArrayRead(da,X,&x);
216:   DMDAVecRestoreArray(da,F,&f);
217:   return(0);
218: }

220: /* --------------------------------------------------------------------- */
221: /*
222:   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
223: */
224: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
225: {
226:   User           user = (User)ptr;
228:   DMDALocalInfo  info;
229:   PetscInt       i;
230:   PetscReal      hx;
231:   DM             da;
232:   Field          *x,*xdot;

235:   TSGetDM(ts,&da);
236:   DMDAGetLocalInfo(da,&info);
237:   hx   = 1.0/(PetscReal)(info.mx-1);

239:   /* Get pointers to vector data */
240:   DMDAVecGetArrayRead(da,X,&x);
241:   DMDAVecGetArrayRead(da,Xdot,&xdot);

243:   /* Compute function over the locally owned part of the grid */
244:   for (i=info.xs; i<info.xs+info.xm; i++) {
245:     if (i == 0 || i == info.mx-1) {
246:       const PetscInt    row        = i,col = i;
247:       const PetscScalar vals[2][2] = {{hx,0},{0,hx}};
248:       MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES);
249:     } else {
250:       const PetscInt    row           = i,col[] = {i-1,i,i+1};
251:       const PetscScalar dxxL          = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
252:       const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
253:                                          {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
254:       MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);
255:     }
256:   }

258:   /* Restore vectors */
259:   DMDAVecRestoreArrayRead(da,X,&x);
260:   DMDAVecRestoreArrayRead(da,Xdot,&xdot);

262:   MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
263:   MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
264:   if (J != Jpre) {
265:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
266:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
267:   }
268:   return(0);
269: }

271: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
272: {
273:   User           user = (User)ctx;
274:   DM             da;
275:   PetscInt       i;
276:   DMDALocalInfo  info;
277:   Field          *x;
278:   PetscReal      hx;

282:   TSGetDM(ts,&da);
283:   DMDAGetLocalInfo(da,&info);
284:   hx   = 1.0/(PetscReal)(info.mx-1);

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

289:   /* Compute function over the locally owned part of the grid */
290:   for (i=info.xs; i<info.xs+info.xm; i++) {
291:     PetscReal xi = i*hx;
292:     x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi);
293:     x[i].v = user->vleft*(1.-xi) + user->vright*xi;
294:   }
295:   DMDAVecRestoreArray(da,X,&x);
296:   return(0);
297: }

299: /*TEST

301:     build:
302:       requires: c99

304:     test:
305:       args: -nox -da_grid_x 20 -ts_monitor_draw_solution -ts_type rosw -ts_rosw_type 2p -ts_dt 5e-2 -ts_adapt_type none

307: TEST*/