Actual source code: ex14.c

  1: /*$Id: ex14.c,v 1.19 2001/04/10 19:37:05 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:    double      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: int main(int argc,char **argv)
 67: {
 68:   SNES                   snes;                 /* nonlinear solver */
 69:   Vec                    x,r;                  /* solution, residual vectors */
 70:   Mat                    J;                    /* Jacobian matrix */
 71:   AppCtx                 user;                 /* user-defined work context */
 72:   int                    its;                  /* iterations for convergence */
 73:   PetscTruth             matrix_free,coloring;
 74:   int                    ierr;
 75:   double                 bratu_lambda_max = 6.81,bratu_lambda_min = 0.,fnorm;
 76:   MatFDColoring          matfdcoloring;

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

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

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

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

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

104:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105:      Extract global vectors from DA; then duplicate for remaining
106:      vectors that are the same types
107:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108:   DACreateGlobalVector(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:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129:   PetscOptionsHasName(PETSC_NULL,"-snes_mf",&matrix_free);
130:   PetscOptionsHasName(PETSC_NULL,"-fdcoloring",&coloring);
131:   if (!matrix_free) {
132:     if (coloring) {
133:       ISColoring    iscoloring;

135:       DAGetColoring(user.da,IS_COLORING_GLOBAL,MATMPIAIJ,&iscoloring,&J);
136:       MatFDColoringCreate(J,iscoloring,&matfdcoloring);
137:       ISColoringDestroy(iscoloring);
138:       MatFDColoringSetFunction(matfdcoloring,(int (*)(void))FormFunction,&user);
139:       MatFDColoringSetFromOptions(matfdcoloring);
140:       SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,matfdcoloring);
141:     } else {
142:       DAGetColoring(user.da,IS_COLORING_GLOBAL,MATMPIAIJ,PETSC_IGNORE,&J);
143:       SNESSetJacobian(snes,J,J,FormJacobian,&user);
144:     }
145:   }

147:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148:      Customize nonlinear solver; set runtime options
149:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150:   SNESSetFromOptions(snes);

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

161:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162:      Solve nonlinear system
163:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164:   SNESSolve(snes,x,&its);

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

173:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174:      Free work space.  All PETSc objects should be destroyed when they
175:      are no longer needed.
176:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

190:   return(0);
191: }
192: /* ------------------------------------------------------------------- */
193: /* 
194:    FormInitialGuess - Forms initial approximation.

196:    Input Parameters:
197:    user - user-defined application context
198:    X - vector

200:    Output Parameter:
201:    X - vector
202:  */
203: int FormInitialGuess(AppCtx *user,Vec X)
204: {
205:   int     i,j,k,Mx,My,Mz,ierr,xs,ys,zs,xm,ym,zm;
206:   double  lambda,temp1,hx,hy,hz,tempk,tempj;
207:   Scalar  ***x;

210:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
211:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

213:   lambda = user->param;
214:   hx     = 1.0/(double)(Mx-1);
215:   hy     = 1.0/(double)(My-1);
216:   hz     = 1.0/(double)(Mz-1);
217:   temp1  = lambda/(lambda + 1.0);

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

228:   /*
229:      Get local grid boundaries (for 3-dimensional DA):
230:        xs, ys, zs   - starting grid indices (no ghost points)
231:        xm, ym, zm   - widths of local grid (no ghost points)

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

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

254:   /*
255:      Restore vector
256:   */
257:   DAVecRestoreArray(user->da,X,(void**)&x);
258:   return(0);
259: }
260: /* ------------------------------------------------------------------- */
261: /* 
262:    FormFunction - Evaluates nonlinear function, F(x).

264:    Input Parameters:
265: .  snes - the SNES context
266: .  X - input vector
267: .  ptr - optional user-defined context, as set by SNESSetFunction()

269:    Output Parameter:
270: .  F - function vector
271:  */
272: int FormFunction(SNES snes,Vec X,Vec F,void *ptr)
273: {
274:   AppCtx  *user = (AppCtx*)ptr;
275:   int     ierr,i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
276:   double  two = 2.0,lambda,hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc;
277:   Scalar  u_north,u_south,u_east,u_west,u_up,u_down,u,u_xx,u_yy,u_zz,***x,***f;
278:   Vec     localX;

281:   DAGetLocalVector(user->da,&localX);
282:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
283:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

285:   lambda = user->param;
286:   hx     = 1.0/(double)(Mx-1);
287:   hy     = 1.0/(double)(My-1);
288:   hz     = 1.0/(double)(Mz-1);
289:   sc     = hx*hy*hz*lambda;
290:   hxhzdhy = hx*hz/hy;
291:   hyhzdhx = hy*hz/hx;
292:   hxhydhz = hx*hy/hz;

294:   /*
295:      Scatter ghost points to local vector,using the 2-step process
296:         DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
297:      By placing code between these two statements, computations can be
298:      done while messages are in transition.
299:   */
300:   DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
301:   DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

303:   /*
304:      Get pointers to vector data
305:   */
306:   DAVecGetArray(user->da,localX,(void**)&x);
307:   DAVecGetArray(user->da,F,(void**)&f);

309:   /*
310:      Get local grid boundaries
311:   */
312:   DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);

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

339:   /*
340:      Restore vectors
341:   */
342:   DAVecRestoreArray(user->da,localX,(void**)&x);
343:   DAVecRestoreArray(user->da,F,(void**)&f);
344:   DARestoreLocalVector(user->da,&localX);
345:   PetscLogFlops(11*ym*xm);
346:   return(0);
347: }
348: /* ------------------------------------------------------------------- */
349: /*
350:    FormJacobian - Evaluates Jacobian matrix.

352:    Input Parameters:
353: .  snes - the SNES context
354: .  x - input vector
355: .  ptr - optional user-defined context, as set by SNESSetJacobian()

357:    Output Parameters:
358: .  A - Jacobian matrix
359: .  B - optionally different preconditioning matrix
360: .  flag - flag indicating matrix structure

362: */
363: int FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
364: {
365:   AppCtx     *user = (AppCtx*)ptr;  /* user-defined application context */
366:   Mat        jac = *B;                /* Jacobian matrix */
367:   Vec        localX;
368:   int        ierr,i,j,k,Mx,My,Mz;
369:   MatStencil col[7],row;
370:   int        xs,ys,zs,xm,ym,zm;
371:   Scalar     lambda,v[7],hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc,***x;


375:   DAGetLocalVector(user->da,&localX);
376:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
377:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

379:   lambda = user->param;
380:   hx     = 1.0/(double)(Mx-1);
381:   hy     = 1.0/(double)(My-1);
382:   hz     = 1.0/(double)(Mz-1);
383:   sc     = hx*hy*hz*lambda;
384:   hxhzdhy = hx*hz/hy;
385:   hyhzdhx = hy*hz/hx;
386:   hxhydhz = hx*hy/hz;

388:   /*
389:      Scatter ghost points to local vector, using the 2-step process
390:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
391:      By placing code between these two statements, computations can be
392:      done while messages are in transition.
393:   */
394:   DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
395:   DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

397:   /*
398:      Get pointer to vector data
399:   */
400:   DAVecGetArray(user->da,localX,(void**)&x);

402:   /*
403:      Get local grid boundaries
404:   */
405:   DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);

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

443:   /* 
444:      Assemble matrix, using the 2-step process:
445:        MatAssemblyBegin(), MatAssemblyEnd().
446:   */
447:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
448:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);

450:   /*
451:      Normally since the matrix has already been assembled above; this
452:      would do nothing. But in the matrix free mode -snes_mf_operator
453:      this tells the "matrix-free" matrix that a new linear system solve
454:      is about to be done.
455:   */

457:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
458:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);

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


480:   /*
481:      Tell the matrix we will never add a new nonzero location to the
482:      matrix. If we do, it will generate an error.
483:   */
484:   MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR);
485:   return(0);
486: }