Actual source code: ex15.c

petsc-3.10.5 2019-03-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 */
 51:   DM             da;
 52:   PetscReal      dt;
 53:   AppCtx         user;              /* user-defined work context */
 54:   SNES           snes;
 55:   PetscInt       Jtype; /* Jacobian type
 56:                             0: user provide Jacobian;
 57:                             1: slow finite difference;
 58:                             2: fd with coloring; */

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

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

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

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

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

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

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

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

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

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

149:   PetscFinalize();
150:   return ierr;
151: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

345: /*TEST

347:     test:
348:       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -ts_monitor

350:     test:
351:       suffix: 2
352:       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 2 -ts_monitor

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

359:     test:
360:       suffix: 4
361:       requires: !single
362:       nsize: 2
363:       args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor 

365:     test:
366:       suffix: 5
367:       nsize: 1
368:       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 1 -ts_monitor 


371: TEST*/