Actual source code: ex15.c

petsc-3.3-p7 2013-05-11
  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


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

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

 18: Program usage:  
 19:    mpiexec -n <procs> ./ex15 [-help] [all PETSc options] 
 20:    e.g., mpiexec -n 2 ./ex15 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor 
 21:          ./ex15 -da_grid_x 40 -da_grid_y 40 -drawcontours -draw_pause .1 -boundary 1 
 22:          ./ex15 -da_grid_x 40 -da_grid_y 40 -drawcontours -draw_pause .1 -boundary 1 -Jtype 2 -nstencilpts 9
 23:         
 24: */

 26: #include <petscdmda.h>
 27: #include <petscts.h>

 29: /* 
 30:    User-defined data structures and routines
 31: */
 32: /* MonitorCtx: used by MyTSMonitor() */
 33: typedef struct {
 34:    PetscBool drawcontours;
 35: } MonitorCtx;

 37: /* AppCtx: used by FormIFunction() and FormIJacobian() */
 38: typedef struct {
 39:   DM             da;
 40:   PetscInt       nstencilpts;    /* number of stencil points: 5 or 9 */
 41:   PetscReal      c;
 42:   PetscInt       boundary;       /* Type of boundary condition */
 43:   PetscBool      viewJacobian;
 44: } AppCtx;

 46: extern PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
 47: extern PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
 48: extern PetscErrorCode FormInitialSolution(Vec,void*);
 49: extern PetscErrorCode MyTSMonitor(TS,PetscInt,PetscReal,Vec,void*);

 53: int main(int argc,char **argv)
 54: {
 55:   TS             ts;                   /* nonlinear solver */
 56:   Vec            u,r;                  /* solution, residual vectors */
 57:   Mat            J,Jmf = PETSC_NULL;   /* Jacobian matrices */
 58:   PetscInt       maxsteps = 1000;      /* iterations for convergence */
 60:   DM             da;
 61:   MatFDColoring  matfdcoloring = PETSC_NULL;
 62:   PetscReal      dt,ftime;
 63:   MonitorCtx     usermonitor;       /* user-defined monitor context */
 64:   AppCtx         user;              /* user-defined work context */
 65:   SNES           snes;
 66:   PetscInt       Jtype; /* Jacobian type
 67:                             0: user provide Jacobian; 
 68:                             1: slow finite difference;
 69:                             2: fd with coloring; */

 71:   PetscInitialize(&argc,&argv,(char *)0,help);
 72:   /* Initialize user application context */
 73:   user.da            = PETSC_NULL;
 74:   user.nstencilpts   = 5;
 75:   user.c             = -30.0;
 76:   user.boundary      = 0; /* 0: Drichlet BC; 1: Neumann BC */
 77:   user.viewJacobian  = PETSC_FALSE;
 78:   PetscOptionsGetInt(PETSC_NULL,"-nstencilpts",&user.nstencilpts,PETSC_NULL);
 79:   PetscOptionsGetInt(PETSC_NULL,"-boundary",&user.boundary,PETSC_NULL);
 80:   PetscOptionsHasName(PETSC_NULL,"-viewJacobian",&user.viewJacobian);

 82:   usermonitor.drawcontours = PETSC_FALSE;
 83:   PetscOptionsHasName(PETSC_NULL,"-drawcontours",&usermonitor.drawcontours);

 85:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 86:      Create distributed array (DMDA) to manage parallel grid and vectors
 87:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 88:   if (user.nstencilpts == 5){
 89:     DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-11,-11,PETSC_DECIDE,PETSC_DECIDE,1,1,PETSC_NULL,PETSC_NULL,&da);
 90:   } else if (user.nstencilpts == 9){
 91:     DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,-11,-11,PETSC_DECIDE,PETSC_DECIDE,1,1,PETSC_NULL,PETSC_NULL,&da);
 92:   } else {
 93:     SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"nstencilpts %d is not supported",user.nstencilpts);
 94:   }
 95:   user.da = da;

 97:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 98:      Extract global vectors from DMDA; 
 99:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100:   DMCreateGlobalVector(da,&u);
101:   VecDuplicate(u,&r);

103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:      Create timestepping solver context
105:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106:   TSCreate(PETSC_COMM_WORLD,&ts);
107:   TSSetProblemType(ts,TS_NONLINEAR);
108:   TSSetType(ts,TSBEULER);
109:   TSSetIFunction(ts,r,FormIFunction,&user);
110:   TSSetDuration(ts,maxsteps,1.0);
111:   TSMonitorSet(ts,MyTSMonitor,&usermonitor,PETSC_NULL);

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:      Set initial conditions
115:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116:   FormInitialSolution(u,&user);
117:   TSSetSolution(ts,u);
118:   dt   = .01;
119:   TSSetInitialTimeStep(ts,0.0,dt);

121:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122:    Set Jacobian evaluation routine 
123:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124:   DMCreateMatrix(da,MATAIJ,&J);
125:   Jtype = 0;
126:   PetscOptionsGetInt(PETSC_NULL, "-Jtype",&Jtype,PETSC_NULL);
127:   if (Jtype == 0){ /* use user provided Jacobian evaluation routine */
128:     if (user.nstencilpts != 5) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"user Jacobian routine FormIJacobian() does not support nstencilpts=%D",user.nstencilpts);
129:     TSSetIJacobian(ts,J,J,FormIJacobian,&user);
130:   } else { /* use finite difference Jacobian J as preconditioner and '-snes_mf_operator' for Mat*vec */
131:     TSGetSNES(ts,&snes);
132:     MatCreateSNESMF(snes,&Jmf);
133:     if (Jtype == 1){ /* slow finite difference J; */
134:       SNESSetJacobian(snes,Jmf,J,SNESDefaultComputeJacobian,PETSC_NULL);
135:     } else if (Jtype == 2){ /* Use coloring to compute  finite difference J efficiently */
136:       ISColoring iscoloring;
137:       DMCreateColoring(da,IS_COLORING_GLOBAL,MATAIJ,&iscoloring);
138:       MatFDColoringCreate(J,iscoloring,&matfdcoloring);
139:       MatFDColoringSetFromOptions(matfdcoloring);
140:       ISColoringDestroy(&iscoloring);
141:       MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
142:       SNESSetJacobian(snes,Jmf,J,SNESDefaultComputeJacobianColor,matfdcoloring);
143:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Jtype is not supported");
144:   }

146:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147:    Sets various TS parameters from user options 
148:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149:   TSSetFromOptions(ts);

151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152:      Solve nonlinear system
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154:   TSSolve(ts,u,&ftime);

156:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157:      Free work space.
158:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159:   MatDestroy(&J);
160:   MatFDColoringDestroy(&matfdcoloring);
161:   MatDestroy(&Jmf);
162:   VecDestroy(&u);
163:   VecDestroy(&r);
164:   TSDestroy(&ts);
165:   DMDestroy(&da);

167:   PetscFinalize();
168:   return(0);
169: }

171: /* --------------------------------------------------------------------- */
172: /*
173:   FormIFunction = Udot - RHSFunction
174: */
177: PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
178: {
180:   AppCtx         *user=(AppCtx*)ctx;
181:   DM             da = (DM)user->da;
182:   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
183:   PetscReal      hx,hy,sx,sy;
184:   PetscScalar    u,uxx,uyy,**uarray,**f,**udot;
185:   Vec            localU;

188:   DMGetLocalVector(da,&localU);
189:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
190:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

192:   hx     = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
193:   hy     = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
194:   if (user->nstencilpts == 9 && hx != hy)SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"hx must equal hy when nstencilpts = 9 for this example");

196:   /*
197:      Scatter ghost points to local vector,using the 2-step process
198:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
199:      By placing code between these two statements, computations can be
200:      done while messages are in transition.
201:   */
202:   DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
203:   DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);

205:   /* Get pointers to vector data */
206:   DMDAVecGetArray(da,localU,&uarray);
207:   DMDAVecGetArray(da,F,&f);
208:   DMDAVecGetArray(da,Udot,&udot);

210:   /* Get local grid boundaries */
211:   DMDAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

213:   /* Compute function over the locally owned part of the grid */
214:   for (j=ys; j<ys+ym; j++) {
215:     for (i=xs; i<xs+xm; i++) {
216:       /* Boundary conditions */
217:       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
218:         if (user->boundary == 0){ /* Drichlet BC */
219:           f[j][i] = uarray[j][i]; /* F = U */
220:         } else {                  /* Neumann BC */
221:           if (i == 0 && j == 0){              /* SW corner */
222:             f[j][i] = uarray[j][i] - uarray[j+1][i+1];
223:           } else if (i == Mx-1 && j == 0){    /* SE corner */
224:             f[j][i] = uarray[j][i] - uarray[j+1][i-1];
225:           } else if (i == 0 && j == My-1){    /* NW corner */
226:             f[j][i] = uarray[j][i] - uarray[j-1][i+1];
227:           } else if (i == Mx-1 && j == My-1){ /* NE corner */
228:             f[j][i] = uarray[j][i] - uarray[j-1][i-1];
229:           } else if (i == 0){                  /* Left */
230:             f[j][i] = uarray[j][i] - uarray[j][i+1];
231:           } else if (i == Mx-1){               /* Right */
232:             f[j][i] = uarray[j][i] - uarray[j][i-1];
233:           } else if (j == 0) {                 /* Bottom */
234:             f[j][i] = uarray[j][i] - uarray[j+1][i];
235:           } else if (j == My-1){               /* Top */
236:             f[j][i] = uarray[j][i] - uarray[j-1][i];
237:           }
238:         }
239:       } else { /* Interior */
240:         u       = uarray[j][i];
241:         /* 5-point stencil */
242:         uxx     = (-2.0*u + uarray[j][i-1] + uarray[j][i+1]);
243:         uyy     = (-2.0*u + uarray[j-1][i] + uarray[j+1][i]);
244:         if (user->nstencilpts == 9){
245:         /* 9-point stencil: assume hx=hy */
246:           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;
247:           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;
248:         }
249:         f[j][i] = udot[j][i] - (uxx*sx + uyy*sy);
250:       }
251:     }
252:   }

254:   /* Restore vectors */
255:   DMDAVecRestoreArray(da,localU,&uarray);
256:   DMDAVecRestoreArray(da,F,&f);
257:   DMDAVecRestoreArray(da,Udot,&udot);
258:   DMRestoreLocalVector(da,&localU);
259:   PetscLogFlops(11.0*ym*xm);
260:   return(0);
261: }

263: /* --------------------------------------------------------------------- */
264: /*
265:   FormIJacobian() - Compute IJacobian = dF/dU + a dF/dUdot 
266:   This routine is not used with option '-use_coloring'
267: */
270: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat *J,Mat *Jpre,MatStructure *str,void *ctx)
271: {
273:   PetscInt       i,j,Mx,My,xs,ys,xm,ym,nc;
274:   AppCtx         *user = (AppCtx*)ctx;
275:   DM             da = (DM)user->da;
276:   MatStencil     col[5],row;
277:   PetscScalar    vals[5],hx,hy,sx,sy;

280:   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);
281:   DMDAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

283:   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
284:   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
285: 
286:   for (j=ys; j<ys+ym; j++){
287:     for (i=xs; i<xs+xm; i++){
288:       nc    = 0;
289:       row.j = j; row.i = i;
290:       if (user->boundary == 0 && (i == 0 || i == Mx-1 || j == 0 || j == My-1)) {
291:         col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0;

293:       } else if (user->boundary > 0 && i == 0) {  /* Left Neumann */
294:         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
295:         col[nc].j = j; col[nc].i = i+1; vals[nc++] = -1.0;
296:       } else if (user->boundary > 0 && i == Mx-1){/* Right Neumann */
297:         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
298:         col[nc].j = j; col[nc].i = i-1; vals[nc++] = -1.0;
299:       } else if (user->boundary > 0 && j == 0) {  /* Bottom Neumann */
300:         col[nc].j = j;   col[nc].i = i; vals[nc++] = 1.0;
301:         col[nc].j = j+1; col[nc].i = i; vals[nc++] = -1.0;
302:       } else if (user->boundary > 0 && j == My-1){/* Top Neumann */
303:         col[nc].j = j;   col[nc].i = i;  vals[nc++] = 1.0;
304:         col[nc].j = j-1; col[nc].i = i;  vals[nc++] = -1.0;
305:       } else {   /* Interior */
306:         col[nc].j = j-1; col[nc].i = i;   vals[nc++] = -sy;
307:         col[nc].j = j;   col[nc].i = i-1; vals[nc++] = -sx;
308:         col[nc].j = j;   col[nc].i = i;   vals[nc++] = 2.0*(sx + sy) + a;
309:         col[nc].j = j;   col[nc].i = i+1; vals[nc++] = -sx;
310:         col[nc].j = j+1; col[nc].i = i;   vals[nc++] = -sy;
311:       }
312:       MatSetValuesStencil(*Jpre,1,&row,nc,col,vals,INSERT_VALUES);
313:     }
314:   }
315:   MatAssemblyBegin(*Jpre,MAT_FINAL_ASSEMBLY);
316:   MatAssemblyEnd(*Jpre,MAT_FINAL_ASSEMBLY);
317:   if (*J != *Jpre) {
318:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
319:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
320:   }

322:   if (user->viewJacobian){
323:     PetscPrintf(((PetscObject)*Jpre)->comm,"Jpre:\n");
324:     MatView(*Jpre,PETSC_VIEWER_STDOUT_WORLD);
325:   }
326:   return(0);
327: }

329: /* ------------------------------------------------------------------- */
332: PetscErrorCode FormInitialSolution(Vec U,void* ptr)
333: {
334:   AppCtx         *user=(AppCtx*)ptr;
335:   DM             da=user->da;
336:   PetscReal      c=user->c;
338:   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
339:   PetscScalar    **u;
340:   PetscReal      hx,hy,x,y,r;

343:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
344:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

346:   hx     = 1.0/(PetscReal)(Mx-1);
347:   hy     = 1.0/(PetscReal)(My-1);

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

352:   /* Get local grid boundaries */
353:   DMDAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

355:   /* Compute function over the locally owned part of the grid */
356:   for (j=ys; j<ys+ym; j++) {
357:     y = j*hy;
358:     for (i=xs; i<xs+xm; i++) {
359:       x = i*hx;
360:       r = PetscSqrtScalar((x-.5)*(x-.5) + (y-.5)*(y-.5));
361:       if (r < .125) {
362:         u[j][i] = PetscExpScalar(c*r*r*r);
363:       } else {
364:         u[j][i] = 0.0;
365:       }
366:     }
367:   }

369:   /* Restore vectors */
370:   DMDAVecRestoreArray(da,U,&u);
371:   return(0);
372: }

376: PetscErrorCode MyTSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ptr)
377: {
379:   PetscReal      norm,vmax,vmin;
380:   MPI_Comm       comm;
381:   MonitorCtx     *user = (MonitorCtx*)ptr;

384:   VecNorm(v,NORM_1,&norm);
385:   VecMax(v,PETSC_NULL,&vmax);
386:   VecMin(v,PETSC_NULL,&vmin);
387:   PetscObjectGetComm((PetscObject)ts,&comm);
388:   PetscPrintf(comm,"timestep %D: time %G, solution norm %G, max %G, min %G\n",step,ptime,norm,vmax,vmin);
389:   if (user->drawcontours){
390:     VecView(v,PETSC_VIEWER_DRAW_WORLD);
391:   }
392:   return(0);
393: }