Actual source code: ex15.c

petsc-master 2017-05-28
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*);

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

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

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

 73:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 74:      Create distributed array (DMDA) to manage parallel grid and vectors
 75:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 76:   if (user.nstencilpts == 5) {
 77:     DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
 78:   } else if (user.nstencilpts == 9) {
 79:     DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
 80:   } else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"nstencilpts %d is not supported",user.nstencilpts);
 81:   DMSetFromOptions(da);
 82:   DMSetUp(da);
 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 ierr;
152: }

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

169:   DMGetLocalVector(da,&localU);
170:   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);

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

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

185:   /* Get pointers to vector data */
186:   DMDAVecGetArrayRead(da,localU,&uarray);
187:   DMDAVecGetArray(da,F,&f);
188:   DMDAVecGetArray(da,Udot,&udot);

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

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

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

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

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

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

264:   for (j=ys; j<ys+ym; j++) {
265:     for (i=xs; i<xs+xm; i++) {
266:       nc    = 0;
267:       row.j = j; row.i = i;
268:       if (user->boundary == 0 && (i == 0 || i == Mx-1 || j == 0 || j == My-1)) {
269:         col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0;

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

300:   if (user->viewJacobian) {
301:     PetscPrintf(PetscObjectComm((PetscObject)Jpre),"Jpre:\n");
302:     MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD);
303:   }
304:   return(0);
305: }

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

319:   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);

321:   hx = 1.0/(PetscReal)(Mx-1);
322:   hy = 1.0/(PetscReal)(My-1);

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

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

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

341:   /* Restore vectors */
342:   DMDAVecRestoreArray(da,U,&u);
343:   return(0);
344: }