Actual source code: ex14.c

petsc-master 2016-07-26
Report Typos and Errors
  2: static char help[] = "Bratu nonlinear PDE in 3d.\n\
  3: We solve the  Bratu (SFI - solid fuel ignition) problem in a 3D rectangular\n\
  4: domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\
  5: The command line options include:\n\
  6:   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
  7:      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n";

  9: /*T
 10:    Concepts: SNES^parallel Bratu example
 11:    Concepts: DMDA^using distributed arrays;
 12:    Processors: n
 13: T*/

 15: /* ------------------------------------------------------------------------

 17:     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 18:     the partial differential equation

 20:             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,

 22:     with boundary conditions

 24:              u = 0  for  x = 0, x = 1, y = 0, y = 1, z = 0, z = 1

 26:     A finite difference approximation with the usual 7-point stencil
 27:     is used to discretize the boundary value problem to obtain a nonlinear
 28:     system of equations.


 31:   ------------------------------------------------------------------------- */

 33: /*
 34:    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
 35:    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
 36:    file automatically includes:
 37:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 38:      petscmat.h - matrices
 39:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 40:      petscviewer.h - viewers               petscpc.h  - preconditioners
 41:      petscksp.h   - linear solvers
 42: */
 43: #include <petscdm.h>
 44: #include <petscdmda.h>
 45: #include <petscsnes.h>


 48: /*
 49:    User-defined application context - contains data needed by the
 50:    application-provided call-back routines, FormJacobian() and
 51:    FormFunction().
 52: */
 53: typedef struct {
 54:   PetscReal param;             /* test problem parameter */
 55:   DM        da;                /* distributed array data structure */
 56: } AppCtx;

 58: /*
 59:    User-defined routines
 60: */
 61: extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*),FormInitialGuess(AppCtx*,Vec);
 62: extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);

 66: int main(int argc,char **argv)
 67: {
 68:   SNES           snes;                         /* nonlinear solver */
 69:   Vec            x,r;                          /* solution, residual vectors */
 70:   Mat            J = NULL;                            /* Jacobian matrix */
 71:   AppCtx         user;                         /* user-defined work context */
 72:   PetscInt       its;                          /* iterations for convergence */
 73:   MatFDColoring  matfdcoloring = NULL;
 74:   PetscBool      matrix_free = PETSC_FALSE,coloring = PETSC_FALSE, coloring_ds = PETSC_FALSE;
 76:   PetscReal      bratu_lambda_max = 6.81,bratu_lambda_min = 0.,fnorm;

 78:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 79:      Initialize program
 80:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 82:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:      Initialize problem parameters
 86:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 87:   user.param = 6.0;
 88:   PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL);
 89:   if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) SETERRQ(PETSC_COMM_SELF,1,"Lambda is out of range");

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 92:      Create nonlinear solver context
 93:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 94:   SNESCreate(PETSC_COMM_WORLD,&snes);

 96:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97:      Create distributed array (DMDA) to manage parallel grid and vectors
 98:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 99:   DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,-4,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,NULL,&user.da);

101:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:      Extract global vectors from DMDA; then duplicate for remaining
103:      vectors that are the same types
104:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105:   DMCreateGlobalVector(user.da,&x);
106:   VecDuplicate(x,&r);

108:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109:      Set function evaluation routine and vector
110:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111:   SNESSetFunction(snes,r,FormFunction,(void*)&user);

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:      Create matrix data structure; set Jacobian evaluation routine

116:      Set Jacobian matrix data structure and default Jacobian evaluation
117:      routine. User can override with:
118:      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
119:                 (unless user explicitly sets preconditioner)
120:      -snes_mf_operator : form preconditioning matrix as set by the user,
121:                          but use matrix-free approx for Jacobian-vector
122:                          products within Newton-Krylov method
123:      -fdcoloring : using finite differences with coloring to compute the Jacobian

125:      Note one can use -matfd_coloring wp or ds the only reason for the -fdcoloring_ds option
126:      below is to test the call to MatFDColoringSetType().
127:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128:   PetscOptionsGetBool(NULL,NULL,"-snes_mf",&matrix_free,NULL);
129:   PetscOptionsGetBool(NULL,NULL,"-fdcoloring",&coloring,NULL);
130:   PetscOptionsGetBool(NULL,NULL,"-fdcoloring_ds",&coloring_ds,NULL);
131:   if (!matrix_free) {
132:     DMSetMatType(user.da,MATAIJ);
133:     DMCreateMatrix(user.da,&J);
134:     if (coloring) {
135:       ISColoring iscoloring;
136:       DMCreateColoring(user.da,IS_COLORING_GLOBAL,&iscoloring);
137:       MatFDColoringCreate(J,iscoloring,&matfdcoloring);
138:       MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormFunction,&user);
139:       if (coloring_ds) {
140:         MatFDColoringSetType(matfdcoloring,MATMFFD_DS);
141:       }
142:       MatFDColoringSetFromOptions(matfdcoloring);
143:       MatFDColoringSetUp(J,iscoloring,matfdcoloring);
144:       SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
145:       ISColoringDestroy(&iscoloring);
146:     } else {
147:       SNESSetJacobian(snes,J,J,FormJacobian,&user);
148:     }
149:   }

151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152:      Customize nonlinear solver; set runtime options
153:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154:   SNESSetDM(snes,user.da);
155:   SNESSetFromOptions(snes);

157:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158:      Evaluate initial guess
159:      Note: The user should initialize the vector, x, with the initial guess
160:      for the nonlinear solver prior to calling SNESSolve().  In particular,
161:      to employ an initial guess of zero, the user should explicitly set
162:      this vector to zero by calling VecSet().
163:   */
164:   FormInitialGuess(&user,x);

166:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167:      Solve nonlinear system
168:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169:   SNESSolve(snes,NULL,x);
170:   SNESGetIterationNumber(snes,&its);

172:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173:      Explicitly check norm of the residual of the solution
174:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175:   FormFunction(snes,x,r,(void*)&user);
176:   VecNorm(r,NORM_2,&fnorm);
177:   PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D fnorm %g\n",its,(double)fnorm);

179:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180:      Free work space.  All PETSc objects should be destroyed when they
181:      are no longer needed.
182:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

184:   MatDestroy(&J);
185:   VecDestroy(&x);
186:   VecDestroy(&r);
187:   SNESDestroy(&snes);
188:   DMDestroy(&user.da);
189:   MatFDColoringDestroy(&matfdcoloring);
190:   PetscFinalize();
191:   return ierr;
192: }
193: /* ------------------------------------------------------------------- */
196: /*
197:    FormInitialGuess - Forms initial approximation.

199:    Input Parameters:
200:    user - user-defined application context
201:    X - vector

203:    Output Parameter:
204:    X - vector
205:  */
206: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
207: {
208:   PetscInt       i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
210:   PetscReal      lambda,temp1,hx,hy,hz,tempk,tempj;
211:   PetscScalar    ***x;

214:   DMDAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

216:   lambda = user->param;
217:   hx     = 1.0/(PetscReal)(Mx-1);
218:   hy     = 1.0/(PetscReal)(My-1);
219:   hz     = 1.0/(PetscReal)(Mz-1);
220:   temp1  = lambda/(lambda + 1.0);

222:   /*
223:      Get a pointer to vector data.
224:        - For default PETSc vectors, VecGetArray() returns a pointer to
225:          the data array.  Otherwise, the routine is implementation dependent.
226:        - You MUST call VecRestoreArray() when you no longer need access to
227:          the array.
228:   */
229:   DMDAVecGetArray(user->da,X,&x);

231:   /*
232:      Get local grid boundaries (for 3-dimensional DMDA):
233:        xs, ys, zs   - starting grid indices (no ghost points)
234:        xm, ym, zm   - widths of local grid (no ghost points)

236:   */
237:   DMDAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);

239:   /*
240:      Compute initial guess over the locally owned part of the grid
241:   */
242:   for (k=zs; k<zs+zm; k++) {
243:     tempk = (PetscReal)(PetscMin(k,Mz-k-1))*hz;
244:     for (j=ys; j<ys+ym; j++) {
245:       tempj = PetscMin((PetscReal)(PetscMin(j,My-j-1))*hy,tempk);
246:       for (i=xs; i<xs+xm; i++) {
247:         if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
248:           /* boundary conditions are all zero Dirichlet */
249:           x[k][j][i] = 0.0;
250:         } else {
251:           x[k][j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,tempj));
252:         }
253:       }
254:     }
255:   }

257:   /*
258:      Restore vector
259:   */
260:   DMDAVecRestoreArray(user->da,X,&x);
261:   return(0);
262: }
263: /* ------------------------------------------------------------------- */
266: /*
267:    FormFunction - Evaluates nonlinear function, F(x).

269:    Input Parameters:
270: .  snes - the SNES context
271: .  X - input vector
272: .  ptr - optional user-defined context, as set by SNESSetFunction()

274:    Output Parameter:
275: .  F - function vector
276:  */
277: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
278: {
279:   AppCtx         *user = (AppCtx*)ptr;
281:   PetscInt       i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
282:   PetscReal      two = 2.0,lambda,hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc;
283:   PetscScalar    u_north,u_south,u_east,u_west,u_up,u_down,u,u_xx,u_yy,u_zz,***x,***f;
284:   Vec            localX;
285:   DM             da;

288:   SNESGetDM(snes,&da);
289:   DMGetLocalVector(da,&localX);
290:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

292:   lambda  = user->param;
293:   hx      = 1.0/(PetscReal)(Mx-1);
294:   hy      = 1.0/(PetscReal)(My-1);
295:   hz      = 1.0/(PetscReal)(Mz-1);
296:   sc      = hx*hy*hz*lambda;
297:   hxhzdhy = hx*hz/hy;
298:   hyhzdhx = hy*hz/hx;
299:   hxhydhz = hx*hy/hz;

301:   /*
302:      Scatter ghost points to local vector,using the 2-step process
303:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
304:      By placing code between these two statements, computations can be
305:      done while messages are in transition.
306:   */
307:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
308:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

310:   /*
311:      Get pointers to vector data
312:   */
313:   DMDAVecGetArrayRead(da,localX,&x);
314:   DMDAVecGetArray(da,F,&f);

316:   /*
317:      Get local grid boundaries
318:   */
319:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);

321:   /*
322:      Compute function over the locally owned part of the grid
323:   */
324:   for (k=zs; k<zs+zm; k++) {
325:     for (j=ys; j<ys+ym; j++) {
326:       for (i=xs; i<xs+xm; i++) {
327:         if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
328:           f[k][j][i] = x[k][j][i];
329:         } else {
330:           u          = x[k][j][i];
331:           u_east     = x[k][j][i+1];
332:           u_west     = x[k][j][i-1];
333:           u_north    = x[k][j+1][i];
334:           u_south    = x[k][j-1][i];
335:           u_up       = x[k+1][j][i];
336:           u_down     = x[k-1][j][i];
337:           u_xx       = (-u_east + two*u - u_west)*hyhzdhx;
338:           u_yy       = (-u_north + two*u - u_south)*hxhzdhy;
339:           u_zz       = (-u_up + two*u - u_down)*hxhydhz;
340:           f[k][j][i] = u_xx + u_yy + u_zz - sc*PetscExpScalar(u);
341:         }
342:       }
343:     }
344:   }

346:   /*
347:      Restore vectors
348:   */
349:   DMDAVecRestoreArrayRead(da,localX,&x);
350:   DMDAVecRestoreArray(da,F,&f);
351:   DMRestoreLocalVector(da,&localX);
352:   PetscLogFlops(11.0*ym*xm);
353:   return(0);
354: }
355: /* ------------------------------------------------------------------- */
358: /*
359:    FormJacobian - Evaluates Jacobian matrix.

361:    Input Parameters:
362: .  snes - the SNES context
363: .  x - input vector
364: .  ptr - optional user-defined context, as set by SNESSetJacobian()

366:    Output Parameters:
367: .  A - Jacobian matrix
368: .  B - optionally different preconditioning matrix

370: */
371: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr)
372: {
373:   AppCtx         *user = (AppCtx*)ptr;  /* user-defined application context */
374:   Vec            localX;
376:   PetscInt       i,j,k,Mx,My,Mz;
377:   MatStencil     col[7],row;
378:   PetscInt       xs,ys,zs,xm,ym,zm;
379:   PetscScalar    lambda,v[7],hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc,***x;
380:   DM             da;

383:   SNESGetDM(snes,&da);
384:   DMGetLocalVector(da,&localX);
385:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

387:   lambda  = user->param;
388:   hx      = 1.0/(PetscReal)(Mx-1);
389:   hy      = 1.0/(PetscReal)(My-1);
390:   hz      = 1.0/(PetscReal)(Mz-1);
391:   sc      = hx*hy*hz*lambda;
392:   hxhzdhy = hx*hz/hy;
393:   hyhzdhx = hy*hz/hx;
394:   hxhydhz = hx*hy/hz;

396:   /*
397:      Scatter ghost points to local vector, using the 2-step process
398:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
399:      By placing code between these two statements, computations can be
400:      done while messages are in transition.
401:   */
402:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
403:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

405:   /*
406:      Get pointer to vector data
407:   */
408:   DMDAVecGetArrayRead(da,localX,&x);

410:   /*
411:      Get local grid boundaries
412:   */
413:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);

415:   /*
416:      Compute entries for the locally owned part of the Jacobian.
417:       - Currently, all PETSc parallel matrix formats are partitioned by
418:         contiguous chunks of rows across the processors.
419:       - Each processor needs to insert only elements that it owns
420:         locally (but any non-local elements will be sent to the
421:         appropriate processor during matrix assembly).
422:       - Here, we set all entries for a particular row at once.
423:       - We can set matrix entries either using either
424:         MatSetValuesLocal() or MatSetValues(), as discussed above.
425:   */
426:   for (k=zs; k<zs+zm; k++) {
427:     for (j=ys; j<ys+ym; j++) {
428:       for (i=xs; i<xs+xm; i++) {
429:         row.k = k; row.j = j; row.i = i;
430:         /* boundary points */
431:         if (i == 0 || j == 0 || k == 0|| i == Mx-1 || j == My-1 || k == Mz-1) {
432:           v[0] = 1.0;
433:           MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
434:         } else {
435:           /* interior grid points */
436:           v[0] = -hxhydhz; col[0].k=k-1;col[0].j=j;  col[0].i = i;
437:           v[1] = -hxhzdhy; col[1].k=k;  col[1].j=j-1;col[1].i = i;
438:           v[2] = -hyhzdhx; col[2].k=k;  col[2].j=j;  col[2].i = i-1;
439:           v[3] = 2.0*(hyhzdhx+hxhzdhy+hxhydhz)-sc*PetscExpScalar(x[k][j][i]);col[3].k=row.k;col[3].j=row.j;col[3].i = row.i;
440:           v[4] = -hyhzdhx; col[4].k=k;  col[4].j=j;  col[4].i = i+1;
441:           v[5] = -hxhzdhy; col[5].k=k;  col[5].j=j+1;col[5].i = i;
442:           v[6] = -hxhydhz; col[6].k=k+1;col[6].j=j;  col[6].i = i;
443:           MatSetValuesStencil(jac,1,&row,7,col,v,INSERT_VALUES);
444:         }
445:       }
446:     }
447:   }
448:   DMDAVecRestoreArrayRead(da,localX,&x);
449:   DMRestoreLocalVector(da,&localX);

451:   /*
452:      Assemble matrix, using the 2-step process:
453:        MatAssemblyBegin(), MatAssemblyEnd().
454:   */
455:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
456:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);

458:   /*
459:      Normally since the matrix has already been assembled above; this
460:      would do nothing. But in the matrix free mode -snes_mf_operator
461:      this tells the "matrix-free" matrix that a new linear system solve
462:      is about to be done.
463:   */

465:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
466:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);

468:   /*
469:      Tell the matrix we will never add a new nonzero location to the
470:      matrix. If we do, it will generate an error.
471:   */
472:   MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
473:   return(0);
474: }