Actual source code: ex14.c

  1: /*$Id: ex14.c,v 1.23 2001/08/07 21:31:17 bsmith Exp $*/

  3: /* Program usage:  mpirun -np <procs> ex14 [-help] [all PETSc options] */

  5: static char help[] = "Bratu nonlinear PDE in 3d.\n\
  6: We solve the  Bratu (SFI - solid fuel ignition) problem in a 3D rectangular\n\
  7: domain, using distributed arrays (DAs) to partition the parallel grid.\n\
  8: The command line options include:\n\
  9:   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
 10:      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n";

 12: /*T
 13:    Concepts: SNES^parallel Bratu example
 14:    Concepts: DA^using distributed arrays;
 15:    Processors: n
 16: T*/

 18: /* ------------------------------------------------------------------------

 20:     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 21:     the partial differential equation
 22:   
 23:             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
 24:   
 25:     with boundary conditions
 26:    
 27:              u = 0  for  x = 0, x = 1, y = 0, y = 1, z = 0, z = 1
 28:   
 29:     A finite difference approximation with the usual 7-point stencil
 30:     is used to discretize the boundary value problem to obtain a nonlinear 
 31:     system of equations.


 34:   ------------------------------------------------------------------------- */

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


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

 60: /* 
 61:    User-defined routines
 62: */
 63: extern int FormFunction(SNES,Vec,Vec,void*),FormInitialGuess(AppCtx*,Vec);
 64: extern int FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);

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

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

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

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

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

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:      Create distributed array (DA) to manage parallel grid and vectors
102:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103:   DACreate3d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,4,4,4,PETSC_DECIDE,PETSC_DECIDE,
104:                     PETSC_DECIDE,1,1,PETSC_NULL,PETSC_NULL,PETSC_NULL,&user.da);

106:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107:      Extract global vectors from DA; then duplicate for remaining
108:      vectors that are the same types
109:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110:   DACreateGlobalVector(user.da,&x);
111:   VecDuplicate(x,&r);

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

118:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119:      Create matrix data structure; set Jacobian evaluation routine

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

130:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131:   PetscOptionsHasName(PETSC_NULL,"-snes_mf",&matrix_free);
132:   PetscOptionsHasName(PETSC_NULL,"-fdcoloring",&coloring);
133:   if (!matrix_free) {
134:     if (coloring) {
135:       ISColoring    iscoloring;

137:       DAGetColoring(user.da,IS_COLORING_LOCAL,&iscoloring);
138:       DAGetMatrix(user.da,MATAIJ,&J);
139:       MatFDColoringCreate(J,iscoloring,&matfdcoloring);
140:       ISColoringDestroy(iscoloring);
141:       MatFDColoringSetFunction(matfdcoloring,(int (*)(void))FormFunction,&user);
142:       MatFDColoringSetFromOptions(matfdcoloring);
143:       SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,matfdcoloring);
144:     } else {
145:       DAGetMatrix(user.da,MATAIJ,&J);
146:       SNESSetJacobian(snes,J,J,FormJacobian,&user);
147:     }
148:   }

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151:      Customize nonlinear solver; set runtime options
152:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153:   SNESSetFromOptions(snes);

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

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

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

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

182:   if (!matrix_free) {
183:     MatDestroy(J);
184:   }
185:   if (coloring) {
186:     MatFDColoringDestroy(matfdcoloring);
187:   }
188:   VecDestroy(x);
189:   VecDestroy(r);
190:   SNESDestroy(snes);
191:   DADestroy(user.da);
192:   PetscFinalize();

194:   return(0);
195: }
196: /* ------------------------------------------------------------------- */
199: /* 
200:    FormInitialGuess - Forms initial approximation.

202:    Input Parameters:
203:    user - user-defined application context
204:    X - vector

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

216:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
217:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

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

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

234:   /*
235:      Get local grid boundaries (for 3-dimensional DA):
236:        xs, ys, zs   - starting grid indices (no ghost points)
237:        xm, ym, zm   - widths of local grid (no ghost points)

239:   */
240:   DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);

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

260:   /*
261:      Restore vector
262:   */
263:   DAVecRestoreArray(user->da,X,&x);
264:   return(0);
265: }
266: /* ------------------------------------------------------------------- */
269: /* 
270:    FormFunction - Evaluates nonlinear function, F(x).

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

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

289:   DAGetLocalVector(user->da,&localX);
290:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
291:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

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

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

311:   /*
312:      Get pointers to vector data
313:   */
314:   DAVecGetArray(user->da,localX,&x);
315:   DAVecGetArray(user->da,F,&f);

317:   /*
318:      Get local grid boundaries
319:   */
320:   DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);

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

347:   /*
348:      Restore vectors
349:   */
350:   DAVecRestoreArray(user->da,localX,&x);
351:   DAVecRestoreArray(user->da,F,&f);
352:   DARestoreLocalVector(user->da,&localX);
353:   PetscLogFlops(11*ym*xm);
354:   return(0);
355: }
356: /* ------------------------------------------------------------------- */
359: /*
360:    FormJacobian - Evaluates Jacobian matrix.

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

367:    Output Parameters:
368: .  A - Jacobian matrix
369: .  B - optionally different preconditioning matrix
370: .  flag - flag indicating matrix structure

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


385:   DAGetLocalVector(user->da,&localX);
386:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
387:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

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

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

407:   /*
408:      Get pointer to vector data
409:   */
410:   DAVecGetArray(user->da,localX,&x);

412:   /*
413:      Get local grid boundaries
414:   */
415:   DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);

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

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

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

467:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
468:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);

470:   /*
471:      Set flag to indicate that the Jacobian matrix retains an identical
472:      nonzero structure throughout all nonlinear iterations (although the
473:      values of the entries change). Thus, we can save some work in setting
474:      up the preconditioner (e.g., no need to redo symbolic factorization for
475:      ILU/ICC preconditioners).
476:       - If the nonzero structure of the matrix is different during
477:         successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
478:         must be used instead.  If you are unsure whether the matrix
479:         structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
480:       - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
481:         believes your assertion and does not check the structure
482:         of the matrix.  If you erroneously claim that the structure
483:         is the same when it actually is not, the new preconditioner
484:         will not function correctly.  Thus, use this optimization
485:         feature with caution!
486:   */
487:   *flag = SAME_NONZERO_PATTERN;


490:   /*
491:      Tell the matrix we will never add a new nonzero location to the
492:      matrix. If we do, it will generate an error.
493:   */
494:   MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR);
495:   return(0);
496: }