Actual source code: ex15.c


  2: static char help[] = "Time-dependent PDE in 2d. Modified from ex13.c for illustrating how to solve DAEs. \n";
  3: /*
  4:    u_t = uxx + uyy
  5:    0 < x < 1, 0 < y < 1;
  6:    At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125
  7:            u(x,y) = 0.0           if r >= .125

  9:    Boundary conditions:
 10:    Drichlet BC:
 11:    At x=0, x=1, y=0, y=1: u = 0.0

 13:    Neumann BC:
 14:    At x=0, x=1: du(x,y,t)/dx = 0
 15:    At y=0, y=1: du(x,y,t)/dy = 0

 17:    mpiexec -n 2 ./ex15 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
 18:          ./ex15 -da_grid_x 40 -da_grid_y 40  -draw_pause .1 -boundary 1 -ts_monitor_draw_solution
 19:          ./ex15 -da_grid_x 40 -da_grid_y 40  -draw_pause .1 -boundary 1 -Jtype 2 -nstencilpts 9

 21: */

 23: #include <petscdm.h>
 24: #include <petscdmda.h>
 25: #include <petscts.h>

 27: /*
 28:    User-defined data structures and routines
 29: */

 31: /* AppCtx: used by FormIFunction() and FormIJacobian() */
 32: typedef struct {
 33:   DM        da;
 34:   PetscInt  nstencilpts;         /* number of stencil points: 5 or 9 */
 35:   PetscReal c;
 36:   PetscInt  boundary;            /* Type of boundary condition */
 37:   PetscBool viewJacobian;
 38: } AppCtx;

 40: extern PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
 41: extern PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
 42: extern PetscErrorCode FormInitialSolution(Vec,void*);

 44: int main(int argc,char **argv)
 45: {
 46:   TS             ts;                   /* nonlinear solver */
 47:   Vec            u,r;                  /* solution, residual vectors */
 48:   Mat            J,Jmf = NULL;   /* Jacobian matrices */
 49:   DM             da;
 50:   PetscReal      dt;
 51:   AppCtx         user;              /* user-defined work context */
 52:   SNES           snes;
 53:   PetscInt       Jtype; /* Jacobian type
 54:                             0: user provide Jacobian;
 55:                             1: slow finite difference;
 56:                             2: fd with coloring; */

 58:   PetscInitialize(&argc,&argv,(char*)0,help);
 59:   /* Initialize user application context */
 60:   user.da           = NULL;
 61:   user.nstencilpts  = 5;
 62:   user.c            = -30.0;
 63:   user.boundary     = 0;  /* 0: Drichlet BC; 1: Neumann BC */
 64:   user.viewJacobian = PETSC_FALSE;

 66:   PetscOptionsGetInt(NULL,NULL,"-nstencilpts",&user.nstencilpts,NULL);
 67:   PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL);
 68:   PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian);

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:      Create distributed array (DMDA) to manage parallel grid and vectors
 72:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 73:   if (user.nstencilpts == 5) {
 74:     DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
 75:   } else if (user.nstencilpts == 9) {
 76:     DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
 77:   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"nstencilpts %d is not supported",user.nstencilpts);
 78:   DMSetFromOptions(da);
 79:   DMSetUp(da);
 80:   user.da = da;

 82:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 83:      Extract global vectors from DMDA;
 84:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 85:   DMCreateGlobalVector(da,&u);
 86:   VecDuplicate(u,&r);

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:      Create timestepping solver context
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 91:   TSCreate(PETSC_COMM_WORLD,&ts);
 92:   TSSetProblemType(ts,TS_NONLINEAR);
 93:   TSSetType(ts,TSBEULER);
 94:   TSSetDM(ts,da);
 95:   TSSetIFunction(ts,r,FormIFunction,&user);
 96:   TSSetMaxTime(ts,1.0);
 97:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

 99:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100:      Set initial conditions
101:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102:   FormInitialSolution(u,&user);
103:   TSSetSolution(ts,u);
104:   dt   = .01;
105:   TSSetTimeStep(ts,dt);

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:    Set Jacobian evaluation routine
109:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110:   DMSetMatType(da,MATAIJ);
111:   DMCreateMatrix(da,&J);
112:   Jtype = 0;
113:   PetscOptionsGetInt(NULL,NULL, "-Jtype",&Jtype,NULL);
114:   if (Jtype == 0) { /* use user provided Jacobian evaluation routine */
116:     TSSetIJacobian(ts,J,J,FormIJacobian,&user);
117:   } else { /* use finite difference Jacobian J as preconditioner and '-snes_mf_operator' for Mat*vec */
118:     TSGetSNES(ts,&snes);
119:     MatCreateSNESMF(snes,&Jmf);
120:     if (Jtype == 1) { /* slow finite difference J; */
121:       SNESSetJacobian(snes,Jmf,J,SNESComputeJacobianDefault,NULL);
122:     } else if (Jtype == 2) { /* Use coloring to compute  finite difference J efficiently */
123:       SNESSetJacobian(snes,Jmf,J,SNESComputeJacobianDefaultColor,0);
124:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Jtype is not supported");
125:   }

127:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:    Sets various TS parameters from user options
129:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130:   TSSetFromOptions(ts);

132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:      Solve nonlinear system
134:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135:   TSSolve(ts,u);

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:      Free work space.
139:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:   MatDestroy(&J);
141:   MatDestroy(&Jmf);
142:   VecDestroy(&u);
143:   VecDestroy(&r);
144:   TSDestroy(&ts);
145:   DMDestroy(&da);

147:   PetscFinalize();
148:   return 0;
149: }

151: /* --------------------------------------------------------------------- */
152: /*
153:   FormIFunction = Udot - RHSFunction
154: */
155: PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
156: {
157:   AppCtx         *user=(AppCtx*)ctx;
158:   DM             da   = (DM)user->da;
159:   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
160:   PetscReal      hx,hy,sx,sy;
161:   PetscScalar    u,uxx,uyy,**uarray,**f,**udot;
162:   Vec            localU;

165:   DMGetLocalVector(da,&localU);
166:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

168:   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
169:   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);

172:   /*
173:      Scatter ghost points to local vector,using the 2-step process
174:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
175:      By placing code between these two statements, computations can be
176:      done while messages are in transition.
177:   */
178:   DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
179:   DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);

181:   /* Get pointers to vector data */
182:   DMDAVecGetArrayRead(da,localU,&uarray);
183:   DMDAVecGetArray(da,F,&f);
184:   DMDAVecGetArray(da,Udot,&udot);

186:   /* Get local grid boundaries */
187:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

189:   /* Compute function over the locally owned part of the grid */
190:   for (j=ys; j<ys+ym; j++) {
191:     for (i=xs; i<xs+xm; i++) {
192:       /* Boundary conditions */
193:       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
194:         if (user->boundary == 0) { /* Drichlet BC */
195:           f[j][i] = uarray[j][i]; /* F = U */
196:         } else {                  /* Neumann BC */
197:           if (i == 0 && j == 0) {              /* SW corner */
198:             f[j][i] = uarray[j][i] - uarray[j+1][i+1];
199:           } else if (i == Mx-1 && j == 0) {    /* SE corner */
200:             f[j][i] = uarray[j][i] - uarray[j+1][i-1];
201:           } else if (i == 0 && j == My-1) {    /* NW corner */
202:             f[j][i] = uarray[j][i] - uarray[j-1][i+1];
203:           } else if (i == Mx-1 && j == My-1) { /* NE corner */
204:             f[j][i] = uarray[j][i] - uarray[j-1][i-1];
205:           } else if (i == 0) {                  /* Left */
206:             f[j][i] = uarray[j][i] - uarray[j][i+1];
207:           } else if (i == Mx-1) {               /* Right */
208:             f[j][i] = uarray[j][i] - uarray[j][i-1];
209:           } else if (j == 0) {                 /* Bottom */
210:             f[j][i] = uarray[j][i] - uarray[j+1][i];
211:           } else if (j == My-1) {               /* Top */
212:             f[j][i] = uarray[j][i] - uarray[j-1][i];
213:           }
214:         }
215:       } else { /* Interior */
216:         u = uarray[j][i];
217:         /* 5-point stencil */
218:         uxx = (-2.0*u + uarray[j][i-1] + uarray[j][i+1]);
219:         uyy = (-2.0*u + uarray[j-1][i] + uarray[j+1][i]);
220:         if (user->nstencilpts == 9) {
221:           /* 9-point stencil: assume hx=hy */
222:           uxx = 2.0*uxx/3.0 + (0.5*(uarray[j-1][i-1]+uarray[j-1][i+1]+uarray[j+1][i-1]+uarray[j+1][i+1]) - 2.0*u)/6.0;
223:           uyy = 2.0*uyy/3.0 + (0.5*(uarray[j-1][i-1]+uarray[j-1][i+1]+uarray[j+1][i-1]+uarray[j+1][i+1]) - 2.0*u)/6.0;
224:         }
225:         f[j][i] = udot[j][i] - (uxx*sx + uyy*sy);
226:       }
227:     }
228:   }

230:   /* Restore vectors */
231:   DMDAVecRestoreArrayRead(da,localU,&uarray);
232:   DMDAVecRestoreArray(da,F,&f);
233:   DMDAVecRestoreArray(da,Udot,&udot);
234:   DMRestoreLocalVector(da,&localU);
235:   PetscLogFlops(11.0*ym*xm);
236:   return 0;
237: }

239: /* --------------------------------------------------------------------- */
240: /*
241:   FormIJacobian() - Compute IJacobian = dF/dU + a dF/dUdot
242:   This routine is not used with option '-use_coloring'
243: */
244: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx)
245: {
246:   PetscInt       i,j,Mx,My,xs,ys,xm,ym,nc;
247:   AppCtx         *user = (AppCtx*)ctx;
248:   DM             da    = (DM)user->da;
249:   MatStencil     col[5],row;
250:   PetscScalar    vals[5],hx,hy,sx,sy;

253:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
254:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

256:   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
257:   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);

259:   for (j=ys; j<ys+ym; j++) {
260:     for (i=xs; i<xs+xm; i++) {
261:       nc    = 0;
262:       row.j = j; row.i = i;
263:       if (user->boundary == 0 && (i == 0 || i == Mx-1 || j == 0 || j == My-1)) {
264:         col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0;

266:       } else if (user->boundary > 0 && i == 0) {  /* Left Neumann */
267:         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
268:         col[nc].j = j; col[nc].i = i+1; vals[nc++] = -1.0;
269:       } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */
270:         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
271:         col[nc].j = j; col[nc].i = i-1; vals[nc++] = -1.0;
272:       } else if (user->boundary > 0 && j == 0) {  /* Bottom Neumann */
273:         col[nc].j = j;   col[nc].i = i; vals[nc++] = 1.0;
274:         col[nc].j = j+1; col[nc].i = i; vals[nc++] = -1.0;
275:       } else if (user->boundary > 0 && j == My-1) { /* Top Neumann */
276:         col[nc].j = j;   col[nc].i = i;  vals[nc++] = 1.0;
277:         col[nc].j = j-1; col[nc].i = i;  vals[nc++] = -1.0;
278:       } else {   /* Interior */
279:         col[nc].j = j-1; col[nc].i = i;   vals[nc++] = -sy;
280:         col[nc].j = j;   col[nc].i = i-1; vals[nc++] = -sx;
281:         col[nc].j = j;   col[nc].i = i;   vals[nc++] = 2.0*(sx + sy) + a;
282:         col[nc].j = j;   col[nc].i = i+1; vals[nc++] = -sx;
283:         col[nc].j = j+1; col[nc].i = i;   vals[nc++] = -sy;
284:       }
285:       MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES);
286:     }
287:   }
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:   }

295:   if (user->viewJacobian) {
296:     PetscPrintf(PetscObjectComm((PetscObject)Jpre),"Jpre:\n");
297:     MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD);
298:   }
299:   return 0;
300: }

302: /* ------------------------------------------------------------------- */
303: PetscErrorCode FormInitialSolution(Vec U,void *ptr)
304: {
305:   AppCtx         *user=(AppCtx*)ptr;
306:   DM             da   =user->da;
307:   PetscReal      c    =user->c;
308:   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
309:   PetscScalar    **u;
310:   PetscReal      hx,hy,x,y,r;

313:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

315:   hx = 1.0/(PetscReal)(Mx-1);
316:   hy = 1.0/(PetscReal)(My-1);

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

321:   /* Get local grid boundaries */
322:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

324:   /* Compute function over the locally owned part of the grid */
325:   for (j=ys; j<ys+ym; j++) {
326:     y = j*hy;
327:     for (i=xs; i<xs+xm; i++) {
328:       x = i*hx;
329:       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
330:       if (r < .125) u[j][i] = PetscExpReal(c*r*r*r);
331:       else u[j][i] = 0.0;
332:     }
333:   }

335:   /* Restore vectors */
336:   DMDAVecRestoreArray(da,U,&u);
337:   return 0;
338: }

340: /*TEST

342:     test:
343:       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -ts_monitor

345:     test:
346:       suffix: 2
347:       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 2 -ts_monitor

349:     test:
350:       suffix: 3
351:       requires: !single
352:       args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor

354:     test:
355:       suffix: 4
356:       requires: !single
357:       nsize: 2
358:       args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor

360:     test:
361:       suffix: 5
362:       nsize: 1
363:       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 1 -ts_monitor

365: TEST*/