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 rectangularn
  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 nonlinearityn
 10:      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)nn";

 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:      petscsles.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*);

 66: #undef __FUNCT__
 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,MATMPIAIJ,&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,MATMPIAIJ,&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,&its);

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

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

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

193:   return(0);
194: }
195: /* ------------------------------------------------------------------- */
196: #undef __FUNCT__
198: /* 
199:    FormInitialGuess - Forms initial approximation.

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

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

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

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

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

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

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

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

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

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

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

288:   DAGetLocalVector(user->da,&localX);
289:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
290:                    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:         DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
304:      By placing code between these two statements, computations can be
305:      done while messages are in transition.
306:   */
307:   DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
308:   DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

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

316:   /*
317:      Get local grid boundaries
318:   */
319:   DAGetCorners(user->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:   DAVecRestoreArray(user->da,localX,(void**)&x);
350:   DAVecRestoreArray(user->da,F,(void**)&f);
351:   DARestoreLocalVector(user->da,&localX);
352:   PetscLogFlops(11*ym*xm);
353:   return(0);
354: }
355: /* ------------------------------------------------------------------- */
356: #undef __FUNCT__
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
369: .  flag - flag indicating matrix structure

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


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

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

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

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

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

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

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

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

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

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


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