Actual source code: ex14.c

petsc-master 2016-12-06
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 FormFunctionLocal(SNES,Vec,Vec,void*);
 62: extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
 63: extern PetscErrorCode FormInitialGuess(AppCtx*,Vec);
 64: extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);

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

 80:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 81:      Initialize program
 82:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

 98:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 99:      Create distributed array (DMDA) to manage parallel grid and vectors
100:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101:   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);
102:   DMSetFromOptions(user.da);
103:   DMSetUp(user.da);
104:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105:      Extract global vectors from DMDA; then duplicate for remaining
106:      vectors that are the same types
107:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108:   DMCreateGlobalVector(user.da,&x);
109:   VecDuplicate(x,&r);

111:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112:      Set function evaluation routine and vector
113:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114:   SNESSetFunction(snes,r,FormFunction,(void*)&user);

116:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117:      Create matrix data structure; set Jacobian evaluation routine

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

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

162:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163:      Customize nonlinear solver; set runtime options
164:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165:   SNESSetDM(snes,user.da);
166:   SNESSetFromOptions(snes);

168:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169:      Evaluate initial guess
170:      Note: The user should initialize the vector, x, with the initial guess
171:      for the nonlinear solver prior to calling SNESSolve().  In particular,
172:      to employ an initial guess of zero, the user should explicitly set
173:      this vector to zero by calling VecSet().
174:   */
175:   FormInitialGuess(&user,x);

177:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178:      Solve nonlinear system
179:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180:   SNESSolve(snes,NULL,x);
181:   SNESGetIterationNumber(snes,&its);

183:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184:      Explicitly check norm of the residual of the solution
185:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186:   FormFunction(snes,x,r,(void*)&user);
187:   VecNorm(r,NORM_2,&fnorm);
188:   PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D fnorm %g\n",its,(double)fnorm);

190:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191:      Free work space.  All PETSc objects should be destroyed when they
192:      are no longer needed.
193:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

195:   MatDestroy(&J);
196:   VecDestroy(&x);
197:   VecDestroy(&r);
198:   SNESDestroy(&snes);
199:   DMDestroy(&user.da);
200:   MatFDColoringDestroy(&matfdcoloring);
201:   PetscFinalize();
202:   return ierr;
203: }
204: /* ------------------------------------------------------------------- */
207: /*
208:    FormInitialGuess - Forms initial approximation.

210:    Input Parameters:
211:    user - user-defined application context
212:    X - vector

214:    Output Parameter:
215:    X - vector
216:  */
217: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
218: {
219:   PetscInt       i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
221:   PetscReal      lambda,temp1,hx,hy,hz,tempk,tempj;
222:   PetscScalar    ***x;

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

227:   lambda = user->param;
228:   hx     = 1.0/(PetscReal)(Mx-1);
229:   hy     = 1.0/(PetscReal)(My-1);
230:   hz     = 1.0/(PetscReal)(Mz-1);
231:   temp1  = lambda/(lambda + 1.0);

233:   /*
234:      Get a pointer to vector data.
235:        - For default PETSc vectors, VecGetArray() returns a pointer to
236:          the data array.  Otherwise, the routine is implementation dependent.
237:        - You MUST call VecRestoreArray() when you no longer need access to
238:          the array.
239:   */
240:   DMDAVecGetArray(user->da,X,&x);

242:   /*
243:      Get local grid boundaries (for 3-dimensional DMDA):
244:        xs, ys, zs   - starting grid indices (no ghost points)
245:        xm, ym, zm   - widths of local grid (no ghost points)

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

250:   /*
251:      Compute initial guess over the locally owned part of the grid
252:   */
253:   for (k=zs; k<zs+zm; k++) {
254:     tempk = (PetscReal)(PetscMin(k,Mz-k-1))*hz;
255:     for (j=ys; j<ys+ym; j++) {
256:       tempj = PetscMin((PetscReal)(PetscMin(j,My-j-1))*hy,tempk);
257:       for (i=xs; i<xs+xm; i++) {
258:         if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
259:           /* boundary conditions are all zero Dirichlet */
260:           x[k][j][i] = 0.0;
261:         } else {
262:           x[k][j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,tempj));
263:         }
264:       }
265:     }
266:   }

268:   /*
269:      Restore vector
270:   */
271:   DMDAVecRestoreArray(user->da,X,&x);
272:   return(0);
273: }
274: /* ------------------------------------------------------------------- */
277: /*
278:    FormFunctionLocal - Evaluates nonlinear function, F(x) on a ghosted patch

280:    Input Parameters:
281: .  snes - the SNES context
282: .  localX - input vector, this contains the ghosted region needed 
283: .  ptr - optional user-defined context, as set by SNESSetFunction()

285:    Output Parameter:
286: .  F - function vector, this does not contain a ghosted region
287:  */
288: PetscErrorCode FormFunctionLocal(SNES snes,Vec localX,Vec F,void *ptr)
289: {
290:   AppCtx         *user = (AppCtx*)ptr;
292:   PetscInt       i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
293:   PetscReal      two = 2.0,lambda,hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc;
294:   PetscScalar    u_north,u_south,u_east,u_west,u_up,u_down,u,u_xx,u_yy,u_zz,***x,***f;
295:   DM             da;

298:   SNESGetDM(snes,&da);
299:   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);

301:   lambda  = user->param;
302:   hx      = 1.0/(PetscReal)(Mx-1);
303:   hy      = 1.0/(PetscReal)(My-1);
304:   hz      = 1.0/(PetscReal)(Mz-1);
305:   sc      = hx*hy*hz*lambda;
306:   hxhzdhy = hx*hz/hy;
307:   hyhzdhx = hy*hz/hx;
308:   hxhydhz = hx*hy/hz;

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:   PetscLogFlops(11.0*ym*xm);
352:   return(0);
353: }
354: /* ------------------------------------------------------------------- */
357: /*
358:    FormFunction - Evaluates nonlinear function, F(x) on the entire domain

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

365:    Output Parameter:
366: .  F - function vector
367:  */
368: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
369: {
371:   Vec            localX;
372:   DM             da;

375:   SNESGetDM(snes,&da);
376:   DMGetLocalVector(da,&localX);

378:   /*
379:      Scatter ghost points to local vector,using the 2-step process
380:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
381:      By placing code between these two statements, computations can be
382:      done while messages are in transition.
383:   */
384:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
385:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

387:   FormFunctionLocal(snes,localX,F,ptr);
388:   DMRestoreLocalVector(da,&localX);
389:   return(0);
390: }
391: /* ------------------------------------------------------------------- */
394: /*
395:    FormJacobian - Evaluates Jacobian matrix.

397:    Input Parameters:
398: .  snes - the SNES context
399: .  x - input vector
400: .  ptr - optional user-defined context, as set by SNESSetJacobian()

402:    Output Parameters:
403: .  A - Jacobian matrix
404: .  B - optionally different preconditioning matrix

406: */
407: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr)
408: {
409:   AppCtx         *user = (AppCtx*)ptr;  /* user-defined application context */
410:   Vec            localX;
412:   PetscInt       i,j,k,Mx,My,Mz;
413:   MatStencil     col[7],row;
414:   PetscInt       xs,ys,zs,xm,ym,zm;
415:   PetscScalar    lambda,v[7],hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc,***x;
416:   DM             da;

419:   SNESGetDM(snes,&da);
420:   DMGetLocalVector(da,&localX);
421:   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);

423:   lambda  = user->param;
424:   hx      = 1.0/(PetscReal)(Mx-1);
425:   hy      = 1.0/(PetscReal)(My-1);
426:   hz      = 1.0/(PetscReal)(Mz-1);
427:   sc      = hx*hy*hz*lambda;
428:   hxhzdhy = hx*hz/hy;
429:   hyhzdhx = hy*hz/hx;
430:   hxhydhz = hx*hy/hz;

432:   /*
433:      Scatter ghost points to local vector, using the 2-step process
434:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
435:      By placing code between these two statements, computations can be
436:      done while messages are in transition.
437:   */
438:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
439:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

441:   /*
442:      Get pointer to vector data
443:   */
444:   DMDAVecGetArrayRead(da,localX,&x);

446:   /*
447:      Get local grid boundaries
448:   */
449:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);

451:   /*
452:      Compute entries for the locally owned part of the Jacobian.
453:       - Currently, all PETSc parallel matrix formats are partitioned by
454:         contiguous chunks of rows across the processors.
455:       - Each processor needs to insert only elements that it owns
456:         locally (but any non-local elements will be sent to the
457:         appropriate processor during matrix assembly).
458:       - Here, we set all entries for a particular row at once.
459:       - We can set matrix entries either using either
460:         MatSetValuesLocal() or MatSetValues(), as discussed above.
461:   */
462:   for (k=zs; k<zs+zm; k++) {
463:     for (j=ys; j<ys+ym; j++) {
464:       for (i=xs; i<xs+xm; i++) {
465:         row.k = k; row.j = j; row.i = i;
466:         /* boundary points */
467:         if (i == 0 || j == 0 || k == 0|| i == Mx-1 || j == My-1 || k == Mz-1) {
468:           v[0] = 1.0;
469:           MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
470:         } else {
471:           /* interior grid points */
472:           v[0] = -hxhydhz; col[0].k=k-1;col[0].j=j;  col[0].i = i;
473:           v[1] = -hxhzdhy; col[1].k=k;  col[1].j=j-1;col[1].i = i;
474:           v[2] = -hyhzdhx; col[2].k=k;  col[2].j=j;  col[2].i = i-1;
475:           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;
476:           v[4] = -hyhzdhx; col[4].k=k;  col[4].j=j;  col[4].i = i+1;
477:           v[5] = -hxhzdhy; col[5].k=k;  col[5].j=j+1;col[5].i = i;
478:           v[6] = -hxhydhz; col[6].k=k+1;col[6].j=j;  col[6].i = i;
479:           MatSetValuesStencil(jac,1,&row,7,col,v,INSERT_VALUES);
480:         }
481:       }
482:     }
483:   }
484:   DMDAVecRestoreArrayRead(da,localX,&x);
485:   DMRestoreLocalVector(da,&localX);

487:   /*
488:      Assemble matrix, using the 2-step process:
489:        MatAssemblyBegin(), MatAssemblyEnd().
490:   */
491:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
492:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);

494:   /*
495:      Normally since the matrix has already been assembled above; this
496:      would do nothing. But in the matrix free mode -snes_mf_operator
497:      this tells the "matrix-free" matrix that a new linear system solve
498:      is about to be done.
499:   */

501:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
502:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);

504:   /*
505:      Tell the matrix we will never add a new nonzero location to the
506:      matrix. If we do, it will generate an error.
507:   */
508:   MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
509:   return(0);
510: }