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 nonlinearity\n\
  8:   -kappa  <parameter>, where <parameter> indicates the problem's nonlinearity\n\
  9:   -mx <xg>, where <xg> = number of grid points in the x-direction\n\
 10:   -my <yg>, where <yg> = number of grid points in the y-direction\n\
 11:   -Nx <npx>, where <npx> = number of processors in the x-direction\n\
 12:   -Ny <npy>, where <npy> = number of processors in the y-direction\n\n";

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

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

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

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

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

111:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112:      Create nonlinear solver context
113:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

115:   SNESCreate(PETSC_COMM_WORLD,&snes);

117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118:      Create vector data structures; set function evaluation routine
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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


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

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

152:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153:      Create matrix data structure; set Jacobian evaluation routine
154:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

188:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189:      Customize nonlinear solver; set runtime options
190:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

211:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212:      Free work space.  All PETSc objects should be destroyed when they
213:      are no longer needed.
214:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

224:   return 0;
225: }
226: /* ------------------------------------------------------------------- */
229: /* 
230:    FormInitialGuess - Forms initial approximation.

232:    Input Parameters:
233:    user - user-defined application context
234:    X - vector

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

246:   mx = user->mx;            my = user->my;            lambda = user->param;
247:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
248:   sc = hx*hy*lambda;        hxdhy = hx/hy;            hydhx = hy/hx;
249:   temp1 = lambda/(lambda + one);

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

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

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

285:   /*
286:      Restore vector
287:   */
288:   VecRestoreArray(localX,&x);

290:   /*
291:      Insert values into global vector
292:   */
293:   DALocalToGlobal(user->da,localX,INSERT_VALUES,X);
294:   return 0;
295: }
296: /* ------------------------------------------------------------------- */
299: /* 
300:    FormFunction - Evaluates nonlinear function, F(x).

302:    Input Parameters:
303: .  snes - the SNES context
304: .  X - input vector
305: .  ptr - optional user-defined context, as set by SNESSetFunction()

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

319:   mx = user->mx;            my = user->my;            lambda = user->param;
320:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
321:   sc = hx*hy*lambda;        hxdhy = hx/hy;            hydhx = hy/hx;
322:   kappa = user->param2;

324:   /*
325:      Scatter ghost points to local vector, using the 2-step process
326:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
327:      By placing code between these two statements, computations can be
328:      done while messages are in transition.
329:   */
330:   DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
331:   DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

333:   /*
334:      Get pointers to vector data
335:   */
336:   VecGetArray(localX,&x);
337:   VecGetArray(localF,&f);

339:   /*
340:      Get local grid boundaries
341:   */
342:   DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
343:   DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL);

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

364:   /*
365:      Restore vectors
366:   */
367:   VecRestoreArray(localX,&x);
368:   VecRestoreArray(localF,&f);

370:   /*
371:      Insert values into global vector
372:   */
373:   DALocalToGlobal(user->da,localF,INSERT_VALUES,F);
374:   PetscLogFlops(11*ym*xm);
375:   return 0;
376: }
377: /* ------------------------------------------------------------------- */
380: /*
381:    FormJacobian - Evaluates Jacobian matrix.

383:    Input Parameters:
384: .  snes - the SNES context
385: .  x - input vector
386: .  ptr - optional user-defined context, as set by SNESSetJacobian()

388:    Output Parameters:
389: .  A - Jacobian matrix
390: .  B - optionally different preconditioning matrix
391: .  flag - flag indicating matrix structure

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

410:   mx = user->mx;            my = user->my;            lambda = user->param;
411:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
412:   sc = hx*hy;               hxdhy = hx/hy;            hydhx = hy/hx;

414:   /*
415:      Scatter ghost points to local vector,using the 2-step process
416:         DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
417:      By placing code between these two statements, computations can be
418:      done while messages are in transition.
419:   */
420:   DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
421:   DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

423:   /*
424:      Get pointer to vector data
425:   */
426:   VecGetArray(localX,&x);

428:   /*
429:      Get local grid boundaries
430:   */
431:   DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
432:   DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL);

434:   /*
435:      Get the global node numbers for all local nodes, including ghost points
436:   */
437:   DAGetGlobalIndices(user->da,&nloc,&ltog);

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

471:   /* 
472:      Assemble matrix, using the 2-step process:
473:        MatAssemblyBegin(), MatAssemblyEnd().
474:      By placing code between these two statements, computations can be
475:      done while messages are in transition.
476:   */
477:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
478:   VecRestoreArray(localX,&x);
479:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);

481:   /*
482:      Set flag to indicate that the Jacobian matrix retains an identical
483:      nonzero structure throughout all nonlinear iterations (although the
484:      values of the entries change). Thus, we can save some work in setting
485:      up the preconditioner (e.g., no need to redo symbolic factorization for
486:      ILU/ICC preconditioners).
487:       - If the nonzero structure of the matrix is different during
488:         successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
489:         must be used instead.  If you are unsure whether the matrix
490:         structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
491:       - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
492:         believes your assertion and does not check the structure
493:         of the matrix.  If you erroneously claim that the structure
494:         is the same when it actually is not, the new preconditioner
495:         will not function correctly.  Thus, use this optimization
496:         feature with caution!
497:   */
498:   *flag = SAME_NONZERO_PATTERN;
499:   /*
500:       Tell the matrix we will never add a new nonzero location to the
501:     matrix. If we do it will generate an error.
502:   */
503:   /* MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR); */
504:   return 0;
505: }