Actual source code: ex15.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  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:    mpiexec -n 2 ./ex15 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
 19:          ./ex15 -da_grid_x 40 -da_grid_y 40  -draw_pause .1 -boundary 1 -ts_monitor_draw_solution
 20:          ./ex15 -da_grid_x 40 -da_grid_y 40  -draw_pause .1 -boundary 1 -Jtype 2 -nstencilpts 9

 22: */

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

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

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

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

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

 63:   PetscInitialize(&argc,&argv,(char*)0,help);
 64:   /* Initialize user application context */
 65:   user.da           = NULL;
 66:   user.nstencilpts  = 5;
 67:   user.c            = -30.0;
 68:   user.boundary     = 0;  /* 0: Drichlet BC; 1: Neumann BC */
 69:   user.viewJacobian = PETSC_FALSE;

 71:   PetscOptionsGetInt(NULL,NULL,"-nstencilpts",&user.nstencilpts,NULL);
 72:   PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL);
 73:   PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian);

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

 85:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 86:      Extract global vectors from DMDA;
 87:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 88:   DMCreateGlobalVector(da,&u);
 89:   VecDuplicate(u,&r);

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

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103:      Set initial conditions
104:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105:   FormInitialSolution(u,&user);
106:   TSSetSolution(ts,u);
107:   dt   = .01;
108:   TSSetInitialTimeStep(ts,0.0,dt);

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

130:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131:    Sets various TS parameters from user options
132:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133:   TSSetFromOptions(ts);

135:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136:      Solve nonlinear system
137:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138:   TSSolve(ts,u);

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

150:   PetscFinalize();
151:   return(0);
152: }

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

171:   DMGetLocalVector(da,&localU);
172:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
173:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

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

179:   /*
180:      Scatter ghost points to local vector,using the 2-step process
181:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
182:      By placing code between these two statements, computations can be
183:      done while messages are in transition.
184:   */
185:   DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
186:   DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);

188:   /* Get pointers to vector data */
189:   DMDAVecGetArrayRead(da,localU,&uarray);
190:   DMDAVecGetArray(da,F,&f);
191:   DMDAVecGetArray(da,Udot,&udot);

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

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

237:   /* Restore vectors */
238:   DMDAVecRestoreArrayRead(da,localU,&uarray);
239:   DMDAVecRestoreArray(da,F,&f);
240:   DMDAVecRestoreArray(da,Udot,&udot);
241:   DMRestoreLocalVector(da,&localU);
242:   PetscLogFlops(11.0*ym*xm);
243:   return(0);
244: }

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

263:   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);
264:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

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

269:   for (j=ys; j<ys+ym; j++) {
270:     for (i=xs; i<xs+xm; i++) {
271:       nc    = 0;
272:       row.j = j; row.i = i;
273:       if (user->boundary == 0 && (i == 0 || i == Mx-1 || j == 0 || j == My-1)) {
274:         col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0;

276:       } else if (user->boundary > 0 && i == 0) {  /* Left Neumann */
277:         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
278:         col[nc].j = j; col[nc].i = i+1; vals[nc++] = -1.0;
279:       } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */
280:         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
281:         col[nc].j = j; col[nc].i = i-1; vals[nc++] = -1.0;
282:       } else if (user->boundary > 0 && j == 0) {  /* Bottom Neumann */
283:         col[nc].j = j;   col[nc].i = i; vals[nc++] = 1.0;
284:         col[nc].j = j+1; col[nc].i = i; vals[nc++] = -1.0;
285:       } else if (user->boundary > 0 && j == My-1) { /* Top Neumann */
286:         col[nc].j = j;   col[nc].i = i;  vals[nc++] = 1.0;
287:         col[nc].j = j-1; col[nc].i = i;  vals[nc++] = -1.0;
288:       } else {   /* Interior */
289:         col[nc].j = j-1; col[nc].i = i;   vals[nc++] = -sy;
290:         col[nc].j = j;   col[nc].i = i-1; vals[nc++] = -sx;
291:         col[nc].j = j;   col[nc].i = i;   vals[nc++] = 2.0*(sx + sy) + a;
292:         col[nc].j = j;   col[nc].i = i+1; vals[nc++] = -sx;
293:         col[nc].j = j+1; col[nc].i = i;   vals[nc++] = -sy;
294:       }
295:       MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES);
296:     }
297:   }
298:   MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
299:   MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
300:   if (J != Jpre) {
301:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
302:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
303:   }

305:   if (user->viewJacobian) {
306:     PetscPrintf(PetscObjectComm((PetscObject)Jpre),"Jpre:\n");
307:     MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD);
308:   }
309:   return(0);
310: }

312: /* ------------------------------------------------------------------- */
315: PetscErrorCode FormInitialSolution(Vec U,void *ptr)
316: {
317:   AppCtx         *user=(AppCtx*)ptr;
318:   DM             da   =user->da;
319:   PetscReal      c    =user->c;
321:   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
322:   PetscScalar    **u;
323:   PetscReal      hx,hy,x,y,r;

326:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
327:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

329:   hx = 1.0/(PetscReal)(Mx-1);
330:   hy = 1.0/(PetscReal)(My-1);

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

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

338:   /* Compute function over the locally owned part of the grid */
339:   for (j=ys; j<ys+ym; j++) {
340:     y = j*hy;
341:     for (i=xs; i<xs+xm; i++) {
342:       x = i*hx;
343:       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
344:       if (r < .125) u[j][i] = PetscExpReal(c*r*r*r);
345:       else u[j][i] = 0.0;
346:     }
347:   }

349:   /* Restore vectors */
350:   DMDAVecRestoreArray(da,U,&u);
351:   return(0);
352: }