Actual source code: ex5.c

  1: /*$Id: ex5.c,v 1.29 2001/08/07 21:31:12 bsmith Exp $*/

  3: static char help[] = "Solves a nonlinear system in parallel with SNES.n
  4: We solve the modified Bratu problem in a 2D rectangular domain,n
  5: using distributed arrays (DAs) to partition the parallel grid.n
  6: The command line options include:n
  7:   -lambda <parameter>, where <parameter> indicates the problem's nonlinearityn
  8:   -kappa  <parameter>, where <parameter> indicates the problem's nonlinearityn
  9:   -mx <xg>, where <xg> = number of grid points in the x-directionn
 10:   -my <yg>, where <yg> = number of grid points in the y-directionn
 11:   -Nx <npx>, where <npx> = number of processors in the x-directionn
 12:   -Ny <npy>, where <npy> = number of processors in the y-directionnn";

 14: /*T
 15:    Concepts: SNES^solving a system of nonlinear equations (parallel Bratu example);
 16:    Concepts: DA^using distributed arrays;
 17:    Processors: n
 18: T*/

 20: /* ------------------------------------------------------------------------

 22:     Modified Solid Fuel Ignition problem.  This problem is modeled by
 23:     the partial differential equation

 25:         -Laplacian u - kappa*PartDer{u}{x} - lambda*exp(u) = 0,

 27:     where

 29:          0 < x,y < 1,
 30:   
 31:     with boundary conditions
 32:    
 33:              u = 0  for  x = 0, x = 1, y = 0, y = 1.
 34:   
 35:     A finite difference approximation with the usual 5-point stencil
 36:     is used to discretize the boundary value problem to obtain a nonlinear 
 37:     system of equations.

 39:   ------------------------------------------------------------------------- */

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

 54: /* 
 55:    User-defined application context - contains data needed by the 
 56:    application-provided call-back routines, FormJacobian() and
 57:    FormFunction().
 58: */
 59: typedef struct {
 60:    PetscReal   param;          /* test problem parameter */
 61:    PetscReal   param2;         /* test problem parameter */
 62:    int         mx,my;          /* discretization in x, y directions */
 63:    Vec         localX,localF; /* ghosted local vector */
 64:    DA          da;             /* distributed array data structure */
 65:    int         rank;           /* processor rank */
 66: } AppCtx;

 68: /* 
 69:    User-defined routines
 70: */
 71: extern int FormFunction(SNES,Vec,Vec,void*),FormInitialGuess(AppCtx*,Vec);
 72: extern int FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);

 74: int main(int argc,char **argv)
 75: {
 76:   SNES       snes;                /* nonlinear solver */
 77:   Vec        x,r;                /* solution, residual vectors */
 78:   Mat        J;                   /* Jacobian matrix */
 79:   AppCtx     user;                /* user-defined work context */
 80:   int        its;                 /* iterations for convergence */
 81:   int        Nx,Ny;              /* number of preocessors in x- and y- directions */
 82:   PetscTruth matrix_free;         /* flag - 1 indicates matrix-free version */
 83:   int        size;                /* number of processors */
 84:   int        m,N,ierr;
 85:   PetscReal  bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
 86:   PetscReal  bratu_kappa_max = 10000,bratu_kappa_min = 0.;

 88:   PetscInitialize(&argc,&argv,(char *)0,help);
 89:   MPI_Comm_rank(PETSC_COMM_WORLD,&user.rank);

 91:   /*
 92:      Initialize problem parameters
 93:   */
 94:   user.mx = 4; user.my = 4; user.param = 6.0; user.param2 = 0.0;
 95:   PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
 96:   PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
 97:   PetscOptionsGetReal(PETSC_NULL,"-lambda",&user.param,PETSC_NULL);
 98:   if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
 99:     SETERRQ(1,"Lambda is out of range");
100:   }
101:   PetscOptionsGetReal(PETSC_NULL,"-kappa",&user.param2,PETSC_NULL);
102:   if (user.param2 >= bratu_kappa_max || user.param2 < bratu_kappa_min) {
103:     SETERRQ(1,"Kappa is out of range");
104:   }
105:   PetscPrintf(PETSC_COMM_WORLD,"Solving the Bratu problem with lambda=%g, kappa=%gn",user.param,user.param2);

107:   N = user.mx*user.my;

109:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110:      Create nonlinear solver context
111:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

113:   SNESCreate(PETSC_COMM_WORLD,&snes);

115:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116:      Create vector data structures; set function evaluation routine
117:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

119:   /*
120:      Create distributed array (DA) to manage parallel grid and vectors
121:   */
122:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
123:   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
124:   PetscOptionsGetInt(PETSC_NULL,"-Nx",&Nx,PETSC_NULL);
125:   PetscOptionsGetInt(PETSC_NULL,"-Ny",&Ny,PETSC_NULL);
126:   if (Nx*Ny != size && (Nx != PETSC_DECIDE || Ny != PETSC_DECIDE))
127:     SETERRQ(1,"Incompatible number of processors:  Nx * Ny != size");
128:   DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,user.mx,user.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&user.da);

130:   /*
131:      Visualize the distribution of the array across the processors
132:   */
133:   /*  DAView(user.da,PETSC_VIEWER_DRAW_WORLD); */


136:   /*
137:      Extract global and local vectors from DA; then duplicate for remaining
138:      vectors that are the same types
139:   */
140:   DACreateGlobalVector(user.da,&x);
141:   DACreateLocalVector(user.da,&user.localX);
142:   VecDuplicate(x,&r);
143:   VecDuplicate(user.localX,&user.localF);

145:   /* 
146:      Set function evaluation routine and vector
147:   */
148:   SNESSetFunction(snes,r,FormFunction,(void*)&user);

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151:      Create matrix data structure; set Jacobian evaluation routine
152:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

154:   /* 
155:      Set Jacobian matrix data structure and default Jacobian evaluation
156:      routine. User can override with:
157:      -snes_fd : default finite differencing approximation of Jacobian
158:      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
159:                 (unless user explicitly sets preconditioner) 
160:      -snes_mf_operator : form preconditioning matrix as set by the user,
161:                          but use matrix-free approx for Jacobian-vector
162:                          products within Newton-Krylov method

164:      Note:  For the parallel case, vectors and matrices MUST be partitioned
165:      accordingly.  When using distributed arrays (DAs) to create vectors,
166:      the DAs determine the problem partitioning.  We must explicitly
167:      specify the local matrix dimensions upon its creation for compatibility
168:      with the vector distribution.  Thus, the generic MatCreate() routine
169:      is NOT sufficient when working with distributed arrays.

171:      Note: Here we only approximately preallocate storage space for the
172:      Jacobian.  See the users manual for a discussion of better techniques
173:      for preallocating matrix memory.
174:   */
175:   PetscOptionsHasName(PETSC_NULL,"-snes_mf",&matrix_free);
176:   if (!matrix_free) {
177:     if (size == 1) {
178:       MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,PETSC_NULL,&J);
179:     } else {
180:       VecGetLocalSize(x,&m);
181:       MatCreateMPIAIJ(PETSC_COMM_WORLD,m,m,N,N,5,PETSC_NULL,3,PETSC_NULL,&J);
182:     }
183:     SNESSetJacobian(snes,J,J,FormJacobian,&user);
184:   }

186:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187:      Customize nonlinear solver; set runtime options
188:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

190:   /*
191:      Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
192:   */
193:   SNESSetFromOptions(snes);

195:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196:      Evaluate initial guess; then solve nonlinear system
197:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198:   /*
199:      Note: The user should initialize the vector, x, with the initial guess
200:      for the nonlinear solver prior to calling SNESSolve().  In particular,
201:      to employ an initial guess of zero, the user should explicitly set
202:      this vector to zero by calling VecSet().
203:   */
204:   FormInitialGuess(&user,x);
205:   SNESSolve(snes,x,&its);
206:   PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %dn",its);

208:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209:      Free work space.  All PETSc objects should be destroyed when they
210:      are no longer needed.
211:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

213:   if (!matrix_free) {
214:     MatDestroy(J);
215:   }
216:   VecDestroy(user.localX); VecDestroy(x);
217:   VecDestroy(user.localF); VecDestroy(r);
218:   SNESDestroy(snes);  DADestroy(user.da);
219:   PetscFinalize();

221:   return 0;
222: }
223: /* ------------------------------------------------------------------- */
224: /* 
225:    FormInitialGuess - Forms initial approximation.

227:    Input Parameters:
228:    user - user-defined application context
229:    X - vector

231:    Output Parameter:
232:    X - vector
233:  */
234: int FormInitialGuess(AppCtx *user,Vec X)
235: {
236:   int          i,j,row,mx,my,ierr,xs,ys,xm,ym,gxm,gym,gxs,gys;
237:   PetscReal    one = 1.0,lambda,temp1,temp,hx,hy,hxdhy,hydhx,sc;
238:   PetscScalar  *x;
239:   Vec          localX = user->localX;

241:   mx = user->mx;            my = user->my;            lambda = user->param;
242:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
243:   sc = hx*hy*lambda;        hxdhy = hx/hy;            hydhx = hy/hx;
244:   temp1 = lambda/(lambda + one);

246:   /*
247:      Get a pointer to vector data.
248:        - For default PETSc vectors,VecGetArray() returns a pointer to
249:          the data array.  Otherwise, the routine is implementation dependent.
250:        - You MUST call VecRestoreArray() when you no longer need access to
251:          the array.
252:   */
253:   VecGetArray(localX,&x);

255:   /*
256:      Get local grid boundaries (for 2-dimensional DA):
257:        xs, ys   - starting grid indices (no ghost points)
258:        xm, ym   - widths of local grid (no ghost points)
259:        gxs, gys - starting grid indices (including ghost points)
260:        gxm, gym - widths of local grid (including ghost points)
261:   */
262:   DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
263:   DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL);

265:   /*
266:      Compute initial guess over the locally owned part of the grid
267:   */
268:   for (j=ys; j<ys+ym; j++) {
269:     temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
270:     for (i=xs; i<xs+xm; i++) {
271:       row = i - gxs + (j - gys)*gxm;
272:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
273:         x[row] = 0.0;
274:         continue;
275:       }
276:       x[row] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
277:     }
278:   }

280:   /*
281:      Restore vector
282:   */
283:   VecRestoreArray(localX,&x);

285:   /*
286:      Insert values into global vector
287:   */
288:   DALocalToGlobal(user->da,localX,INSERT_VALUES,X);
289:   return 0;
290: }
291: /* ------------------------------------------------------------------- */
292: /* 
293:    FormFunction - Evaluates nonlinear function, F(x).

295:    Input Parameters:
296: .  snes - the SNES context
297: .  X - input vector
298: .  ptr - optional user-defined context, as set by SNESSetFunction()

300:    Output Parameter:
301: .  F - function vector
302:  */
303: int FormFunction(SNES snes,Vec X,Vec F,void *ptr)
304: {
305:   AppCtx      *user = (AppCtx*)ptr;
306:   int         ierr,i,j,row,mx,my,xs,ys,xm,ym,gxs,gys,gxm,gym;
307:   PetscReal   two = 2.0,one = 1.0,half = 0.5;
308:   PetscReal   lambda,hx,hy,hxdhy,hydhx,sc;
309:   PetscScalar u,ux,uxx,uyy,*x,*f,kappa;
310:   Vec         localX = user->localX,localF = user->localF;

312:   mx = user->mx;            my = user->my;            lambda = user->param;
313:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
314:   sc = hx*hy*lambda;        hxdhy = hx/hy;            hydhx = hy/hx;
315:   kappa = user->param2;

317:   /*
318:      Scatter ghost points to local vector, using the 2-step process
319:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
320:      By placing code between these two statements, computations can be
321:      done while messages are in transition.
322:   */
323:   DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
324:   DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

326:   /*
327:      Get pointers to vector data
328:   */
329:   VecGetArray(localX,&x);
330:   VecGetArray(localF,&f);

332:   /*
333:      Get local grid boundaries
334:   */
335:   DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
336:   DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL);

338:   /*
339:      Compute function over the locally owned part of the grid
340:   */
341:   for (j=ys; j<ys+ym; j++) {
342:     row = (j - gys)*gxm + xs - gxs - 1;
343:     for (i=xs; i<xs+xm; i++) {
344:       row++;
345:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
346:         f[row] = x[row];
347:         continue;
348:       }
349:       u = x[row];
350:       ux  = (x[row+1] - x[row-1])*half*hy;
351:       uxx = (two*u - x[row-1] - x[row+1])*hydhx;
352:       uyy = (two*u - x[row-gxm] - x[row+gxm])*hxdhy;
353:       f[row] = uxx + uyy - kappa*ux - sc*exp(u);
354:     }
355:   }

357:   /*
358:      Restore vectors
359:   */
360:   VecRestoreArray(localX,&x);
361:   VecRestoreArray(localF,&f);

363:   /*
364:      Insert values into global vector
365:   */
366:   DALocalToGlobal(user->da,localF,INSERT_VALUES,F);
367:   PetscLogFlops(11*ym*xm);
368:   return 0;
369: }
370: /* ------------------------------------------------------------------- */
371: /*
372:    FormJacobian - Evaluates Jacobian matrix.

374:    Input Parameters:
375: .  snes - the SNES context
376: .  x - input vector
377: .  ptr - optional user-defined context, as set by SNESSetJacobian()

379:    Output Parameters:
380: .  A - Jacobian matrix
381: .  B - optionally different preconditioning matrix
382: .  flag - flag indicating matrix structure

384:    Notes:
385:    Due to grid point reordering with DAs, we must always work
386:    with the local grid points, and then transform them to the new
387:    global numbering with the "ltog" mapping (via DAGetGlobalIndices()).
388:    We cannot work directly with the global numbers for the original
389:    uniprocessor grid!
390: */
391: int FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
392: {
393:   AppCtx  *user = (AppCtx*)ptr;  /* user-defined application context */
394:   Mat     jac = *B;                /* Jacobian matrix */
395:   Vec     localX = user->localX;   /* local vector */
396:   int     *ltog;                   /* local-to-global mapping */
397:   int     ierr,i,j,row,mx,my,col[5];
398:   int     nloc,xs,ys,xm,ym,gxs,gys,gxm,gym,grow;
399:   PetscScalar  two = 2.0,one = 1.0,lambda,v[5],hx,hy,hxdhy,hydhx,sc,*x;

401:   mx = user->mx;            my = user->my;            lambda = user->param;
402:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
403:   sc = hx*hy;               hxdhy = hx/hy;            hydhx = hy/hx;

405:   /*
406:      Scatter ghost points to local vector,using the 2-step process
407:         DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
408:      By placing code between these two statements, computations can be
409:      done while messages are in transition.
410:   */
411:   DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
412:   DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

414:   /*
415:      Get pointer to vector data
416:   */
417:   VecGetArray(localX,&x);

419:   /*
420:      Get local grid boundaries
421:   */
422:   DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
423:   DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL);

425:   /*
426:      Get the global node numbers for all local nodes, including ghost points
427:   */
428:   DAGetGlobalIndices(user->da,&nloc,&ltog);

430:   /* 
431:      Compute entries for the locally owned part of the Jacobian.
432:       - Currently, all PETSc parallel matrix formats are partitioned by
433:         contiguous chunks of rows across the processors. The "grow"
434:         parameter computed below specifies the global row number 
435:         corresponding to each local grid point.
436:       - Each processor needs to insert only elements that it owns
437:         locally (but any non-local elements will be sent to the
438:         appropriate processor during matrix assembly). 
439:       - Always specify global row and columns of matrix entries.
440:       - Here, we set all entries for a particular row at once.
441:   */
442:   for (j=ys; j<ys+ym; j++) {
443:     row = (j - gys)*gxm + xs - gxs - 1;
444:     for (i=xs; i<xs+xm; i++) {
445:       row++;
446:       grow = ltog[row];
447:       /* boundary points */
448:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
449:         MatSetValues(jac,1,&grow,1,&grow,&one,INSERT_VALUES);
450:         continue;
451:       }
452:       /* interior grid points */
453:       v[0] = -hxdhy; col[0] = ltog[row - gxm];
454:       v[1] = -hydhx; col[1] = ltog[row - 1];
455:       v[2] = two*(hydhx + hxdhy) - sc*lambda*exp(x[row]); col[2] = grow;
456:       v[3] = -hydhx; col[3] = ltog[row + 1];
457:       v[4] = -hxdhy; col[4] = ltog[row + gxm];
458:       MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
459:     }
460:   }

462:   /* 
463:      Assemble matrix, using the 2-step process:
464:        MatAssemblyBegin(), MatAssemblyEnd().
465:      By placing code between these two statements, computations can be
466:      done while messages are in transition.
467:   */
468:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
469:   VecRestoreArray(localX,&x);
470:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);

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