Actual source code: ex14.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: static char help[] = "Solves a nonlinear system in parallel with a user-defined Newton method.\n\
  3: Uses KSP to solve the linearized Newton sytems.  This solver\n\
  4: is a very simplistic inexact Newton method.  The intent of this code is to\n\
  5: demonstrate the repeated solution of linear sytems with the same nonzero pattern.\n\
  6: \n\
  7: This is NOT the recommended approach for solving nonlinear problems with PETSc!\n\
  8: We urge users to employ the SNES component for solving nonlinear problems whenever\n\
  9: possible, as it offers many advantages over coding nonlinear solvers independently.\n\
 10: \n\
 11: We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\
 12: domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\
 13: The command line options include:\n\
 14:   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
 15:      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\
 16:   -mx <xg>, where <xg> = number of grid points in the x-direction\n\
 17:   -my <yg>, where <yg> = number of grid points in the y-direction\n\
 18:   -Nx <npx>, where <npx> = number of processors in the x-direction\n\
 19:   -Ny <npy>, where <npy> = number of processors in the y-direction\n\n";

 21: /*T
 22:    Concepts: KSP^writing a user-defined nonlinear solver (parallel Bratu example);
 23:    Concepts: DMDA^using distributed arrays;
 24:    Processors: n
 25: T*/

 27: /* ------------------------------------------------------------------------

 29:     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 30:     the partial differential equation

 32:             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,

 34:     with boundary conditions

 36:              u = 0  for  x = 0, x = 1, y = 0, y = 1.

 38:     A finite difference approximation with the usual 5-point stencil
 39:     is used to discretize the boundary value problem to obtain a nonlinear
 40:     system of equations.

 42:     The SNES version of this problem is:  snes/examples/tutorials/ex5.c
 43:     We urge users to employ the SNES component for solving nonlinear
 44:     problems whenever possible, as it offers many advantages over coding
 45:     nonlinear solvers independently.

 47:   ------------------------------------------------------------------------- */

 49: /*
 50:    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
 51:    Include "petscksp.h" so that we can use KSP solvers.  Note that this
 52:    file automatically includes:
 53:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 54:      petscmat.h - matrices
 55:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 56:      petscviewer.h - viewers               petscpc.h  - preconditioners
 57: */
 58: #include <petscdm.h>
 59: #include <petscdmda.h>
 60: #include <petscksp.h>

 62: /*
 63:    User-defined application context - contains data needed by the
 64:    application-provided call-back routines, ComputeJacobian() and
 65:    ComputeFunction().
 66: */
 67: typedef struct {
 68:   PetscReal param;             /* test problem parameter */
 69:   PetscInt  mx,my;             /* discretization in x,y directions */
 70:   Vec       localX;           /* ghosted local vector */
 71:   DM        da;                /* distributed array data structure */
 72: } AppCtx;

 74: /*
 75:    User-defined routines
 76: */
 77: extern PetscErrorCode ComputeFunction(AppCtx*,Vec,Vec),FormInitialGuess(AppCtx*,Vec);
 78: extern PetscErrorCode ComputeJacobian(AppCtx*,Vec,Mat);

 82: int main(int argc,char **argv)
 83: {
 84:   /* -------------- Data to define application problem ---------------- */
 85:   MPI_Comm       comm;                /* communicator */
 86:   KSP            ksp;                /* linear solver */
 87:   Vec            X,Y,F;             /* solution, update, residual vectors */
 88:   Mat            J;                   /* Jacobian matrix */
 89:   AppCtx         user;                /* user-defined work context */
 90:   PetscInt       Nx,Ny;              /* number of preocessors in x- and y- directions */
 91:   PetscMPIInt    size;                /* number of processors */
 92:   PetscReal      bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
 93:   PetscInt       m,N;

 96:   /* --------------- Data to define nonlinear solver -------------- */
 97:   PetscReal    rtol = 1.e-8;          /* relative convergence tolerance */
 98:   PetscReal    xtol = 1.e-8;          /* step convergence tolerance */
 99:   PetscReal    ttol;                  /* convergence tolerance */
100:   PetscReal    fnorm,ynorm,xnorm;     /* various vector norms */
101:   PetscInt     max_nonlin_its = 10;   /* maximum number of iterations for nonlinear solver */
102:   PetscInt     max_functions  = 50;   /* maximum number of function evaluations */
103:   PetscInt     lin_its;               /* number of linear solver iterations for each step */
104:   PetscInt     i;                     /* nonlinear solve iteration number */
105:   PetscBool    no_output = PETSC_FALSE;             /* flag indicating whether to surpress output */

107:   PetscInitialize(&argc,&argv,(char*)0,help);
108:   comm = PETSC_COMM_WORLD;
109:   PetscOptionsGetBool(NULL,"-no_output",&no_output,NULL);

111:   /*
112:      Initialize problem parameters
113:   */
114:   user.mx = 4; user.my = 4; user.param = 6.0;

116:   PetscOptionsGetInt(NULL,"-mx",&user.mx,NULL);
117:   PetscOptionsGetInt(NULL,"-my",&user.my,NULL);
118:   PetscOptionsGetReal(NULL,"-par",&user.param,NULL);
119:   if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) SETERRQ(PETSC_COMM_WORLD,1,"Lambda is out of range");
120:   N = user.mx*user.my;

122:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123:      Create linear solver context
124:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

126:   KSPCreate(comm,&ksp);

128:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129:      Create vector data structures
130:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

132:   /*
133:      Create distributed array (DMDA) to manage parallel grid and vectors
134:   */
135:   MPI_Comm_size(comm,&size);
136:   Nx   = PETSC_DECIDE; Ny = PETSC_DECIDE;
137:   PetscOptionsGetInt(NULL,"-Nx",&Nx,NULL);
138:   PetscOptionsGetInt(NULL,"-Ny",&Ny,NULL);
139:   if (Nx*Ny != size && (Nx != PETSC_DECIDE || Ny != PETSC_DECIDE)) SETERRQ(PETSC_COMM_WORLD,1,"Incompatible number of processors:  Nx * Ny != size");
140:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,&user.da);

142:   /*
143:      Extract global and local vectors from DMDA; then duplicate for remaining
144:      vectors that are the same types
145:   */
146:   DMCreateGlobalVector(user.da,&X);
147:   DMCreateLocalVector(user.da,&user.localX);
148:   VecDuplicate(X,&F);
149:   VecDuplicate(X,&Y);


152:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153:      Create matrix data structure for Jacobian
154:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155:   /*
156:      Note:  For the parallel case, vectors and matrices MUST be partitioned
157:      accordingly.  When using distributed arrays (DMDAs) to create vectors,
158:      the DMDAs determine the problem partitioning.  We must explicitly
159:      specify the local matrix dimensions upon its creation for compatibility
160:      with the vector distribution.  Thus, the generic MatCreate() routine
161:      is NOT sufficient when working with distributed arrays.

163:      Note: Here we only approximately preallocate storage space for the
164:      Jacobian.  See the users manual for a discussion of better techniques
165:      for preallocating matrix memory.
166:   */
167:   if (size == 1) {
168:     MatCreateSeqAIJ(comm,N,N,5,NULL,&J);
169:   } else {
170:     VecGetLocalSize(X,&m);
171:     MatCreateAIJ(comm,m,m,N,N,5,NULL,3,NULL,&J);
172:   }

174:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175:      Customize linear solver; set runtime options
176:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

178:   /*
179:      Set runtime options (e.g.,-ksp_monitor -ksp_rtol <rtol> -ksp_type <type>)
180:   */
181:   KSPSetFromOptions(ksp);

183:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184:      Evaluate initial guess
185:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

187:   FormInitialGuess(&user,X);
188:   ComputeFunction(&user,X,F);   /* Compute F(X)    */
189:   VecNorm(F,NORM_2,&fnorm);     /* fnorm = || F || */
190:   ttol = fnorm*rtol;
191:   if (!no_output) PetscPrintf(comm,"Initial function norm = %g\n",(double)fnorm);

193:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194:      Solve nonlinear system with a user-defined method
195:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

197:   /*
198:       This solver is a very simplistic inexact Newton method, with no
199:       no damping strategies or bells and whistles. The intent of this code
200:       is  merely to demonstrate the repeated solution with KSP of linear
201:       sytems with the same nonzero structure.

203:       This is NOT the recommended approach for solving nonlinear problems
204:       with PETSc!  We urge users to employ the SNES component for solving
205:       nonlinear problems whenever possible with application codes, as it
206:       offers many advantages over coding nonlinear solvers independently.
207:    */

209:   for (i=0; i<max_nonlin_its; i++) {

211:     /*
212:         Compute the Jacobian matrix.  
213:      */
214:     ComputeJacobian(&user,X,J);

216:     /*
217:         Solve J Y = F, where J is the Jacobian matrix.
218:           - First, set the KSP linear operators.  Here the matrix that
219:             defines the linear system also serves as the preconditioning
220:             matrix.
221:           - Then solve the Newton system.
222:      */
223:     KSPSetOperators(ksp,J,J);
224:     KSPSolve(ksp,F,Y);
225:     KSPGetIterationNumber(ksp,&lin_its);

227:     /*
228:        Compute updated iterate
229:      */
230:     VecNorm(Y,NORM_2,&ynorm);       /* ynorm = || Y || */
231:     VecAYPX(Y,-1.0,X);              /* Y <- X - Y      */
232:     VecCopy(Y,X);                   /* X <- Y          */
233:     VecNorm(X,NORM_2,&xnorm);       /* xnorm = || X || */
234:     if (!no_output) {
235:       PetscPrintf(comm,"   linear solve iterations = %D, xnorm=%g, ynorm=%g\n",lin_its,(double)xnorm,(double)ynorm);
236:     }

238:     /*
239:        Evaluate new nonlinear function
240:      */
241:     ComputeFunction(&user,X,F);     /* Compute F(X)    */
242:     VecNorm(F,NORM_2,&fnorm);       /* fnorm = || F || */
243:     if (!no_output) {
244:       PetscPrintf(comm,"Iteration %D, function norm = %g\n",i+1,(double)fnorm);
245:     }

247:     /*
248:        Test for convergence
249:      */
250:     if (fnorm <= ttol) {
251:       if (!no_output) {
252:         PetscPrintf(comm,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)ttol);
253:       }
254:       break;
255:     }
256:     if (ynorm < xtol*(xnorm)) {
257:       if (!no_output) {
258:         PetscPrintf(comm,"Converged due to small update length: %g < %g * %g\n",(double)ynorm,(double)xtol,(double)xnorm);
259:       }
260:       break;
261:     }
262:     if (i > max_functions) {
263:       if (!no_output) {
264:         PetscPrintf(comm,"Exceeded maximum number of function evaluations: %D > %D\n",i,max_functions);
265:       }
266:       break;
267:     }
268:   }
269:   PetscPrintf(comm,"Number of SNES iterations = %D\n",i+1);

271:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272:      Free work space.  All PETSc objects should be destroyed when they
273:      are no longer needed.
274:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

276:   MatDestroy(&J);           VecDestroy(&Y);
277:   VecDestroy(&user.localX); VecDestroy(&X);
278:   VecDestroy(&F);
279:   KSPDestroy(&ksp);  DMDestroy(&user.da);
280:   PetscFinalize();

282:   return 0;
283: }
284: /* ------------------------------------------------------------------- */
287: /*
288:    FormInitialGuess - Forms initial approximation.

290:    Input Parameters:
291:    user - user-defined application context
292:    X - vector

294:    Output Parameter:
295:    X - vector
296:  */
297: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
298: {
299:   PetscInt    i,j,row,mx,my,ierr,xs,ys,xm,ym,gxm,gym,gxs,gys;
300:   PetscReal   one = 1.0,lambda,temp1,temp,hx,hy;
301:   PetscScalar *x;

303:   mx    = user->mx;            my = user->my;            lambda = user->param;
304:   hx    = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
305:   temp1 = lambda/(lambda + one);

307:   /*
308:      Get a pointer to vector data.
309:        - For default PETSc vectors, VecGetArray() returns a pointer to
310:          the data array.  Otherwise, the routine is implementation dependent.
311:        - You MUST call VecRestoreArray() when you no longer need access to
312:          the array.
313:   */
314:   VecGetArray(X,&x);

316:   /*
317:      Get local grid boundaries (for 2-dimensional DMDA):
318:        xs, ys   - starting grid indices (no ghost points)
319:        xm, ym   - widths of local grid (no ghost points)
320:        gxs, gys - starting grid indices (including ghost points)
321:        gxm, gym - widths of local grid (including ghost points)
322:   */
323:   DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
324:   DMDAGetGhostCorners(user->da,&gxs,&gys,NULL,&gxm,&gym,NULL);

326:   /*
327:      Compute initial guess over the locally owned part of the grid
328:   */
329:   for (j=ys; j<ys+ym; j++) {
330:     temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
331:     for (i=xs; i<xs+xm; i++) {
332:       row = i - gxs + (j - gys)*gxm;
333:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
334:         x[row] = 0.0;
335:         continue;
336:       }
337:       x[row] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
338:     }
339:   }

341:   /*
342:      Restore vector
343:   */
344:   VecRestoreArray(X,&x);
345:   return 0;
346: }
347: /* ------------------------------------------------------------------- */
350: /*
351:    ComputeFunction - Evaluates nonlinear function, F(x).

353:    Input Parameters:
354: .  X - input vector
355: .  user - user-defined application context

357:    Output Parameter:
358: .  F - function vector
359:  */
360: PetscErrorCode ComputeFunction(AppCtx *user,Vec X,Vec F)
361: {
363:   PetscInt       i,j,row,mx,my,xs,ys,xm,ym,gxs,gys,gxm,gym;
364:   PetscReal      two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx,sc;
365:   PetscScalar    u,uxx,uyy,*x,*f;
366:   Vec            localX = user->localX;

368:   mx = user->mx;            my = user->my;            lambda = user->param;
369:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
370:   sc = hx*hy*lambda;        hxdhy = hx/hy;            hydhx = hy/hx;

372:   /*
373:      Scatter ghost points to local vector, using the 2-step process
374:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
375:      By placing code between these two statements, computations can be
376:      done while messages are in transition.
377:   */
378:   DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
379:   DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

381:   /*
382:      Get pointers to vector data
383:   */
384:   VecGetArray(localX,&x);
385:   VecGetArray(F,&f);

387:   /*
388:      Get local grid boundaries
389:   */
390:   DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
391:   DMDAGetGhostCorners(user->da,&gxs,&gys,NULL,&gxm,&gym,NULL);

393:   /*
394:      Compute function over the locally owned part of the grid
395:   */
396:   for (j=ys; j<ys+ym; j++) {
397:     row = (j - gys)*gxm + xs - gxs - 1;
398:     for (i=xs; i<xs+xm; i++) {
399:       row++;
400:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
401:         f[row] = x[row];
402:         continue;
403:       }
404:       u      = x[row];
405:       uxx    = (two*u - x[row-1] - x[row+1])*hydhx;
406:       uyy    = (two*u - x[row-gxm] - x[row+gxm])*hxdhy;
407:       f[row] = uxx + uyy - sc*PetscExpScalar(u);
408:     }
409:   }

411:   /*
412:      Restore vectors
413:   */
414:   VecRestoreArray(localX,&x);
415:   VecRestoreArray(F,&f);
416:   PetscLogFlops(11.0*ym*xm);
417:   return 0;
418: }
419: /* ------------------------------------------------------------------- */
422: /*
423:    ComputeJacobian - Evaluates Jacobian matrix.

425:    Input Parameters:
426: .  x - input vector
427: .  user - user-defined application context

429:    Output Parameters:
430: .  jac - Jacobian matrix
431: .  flag - flag indicating matrix structure

433:    Notes:
434:    Due to grid point reordering with DMDAs, we must always work
435:    with the local grid points, and then transform them to the new
436:    global numbering with the "ltog" mapping 
437:    We cannot work directly with the global numbers for the original
438:    uniprocessor grid!
439: */
440: PetscErrorCode ComputeJacobian(AppCtx *user,Vec X,Mat jac)
441: {
442:   PetscErrorCode         ierr;
443:   Vec                    localX = user->localX;   /* local vector */
444:   const PetscInt         *ltog;                   /* local-to-global mapping */
445:   PetscInt               i,j,row,mx,my,col[5];
446:   PetscInt               xs,ys,xm,ym,gxs,gys,gxm,gym,grow;
447:   PetscScalar            two = 2.0,one = 1.0,lambda,v[5],hx,hy,hxdhy,hydhx,sc,*x;
448:   ISLocalToGlobalMapping ltogm;

450:   mx = user->mx;            my = user->my;            lambda = user->param;
451:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
452:   sc = hx*hy;               hxdhy = hx/hy;            hydhx = hy/hx;

454:   /*
455:      Scatter ghost points to local vector, using the 2-step process
456:         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
457:      By placing code between these two statements, computations can be
458:      done while messages are in transition.
459:   */
460:   DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
461:   DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

463:   /*
464:      Get pointer to vector data
465:   */
466:   VecGetArray(localX,&x);

468:   /*
469:      Get local grid boundaries
470:   */
471:   DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
472:   DMDAGetGhostCorners(user->da,&gxs,&gys,NULL,&gxm,&gym,NULL);

474:   /*
475:      Get the global node numbers for all local nodes, including ghost points
476:   */
477:   DMGetLocalToGlobalMapping(user->da,&ltogm);
478:   ISLocalToGlobalMappingGetIndices(ltogm,&ltog);

480:   /*
481:      Compute entries for the locally owned part of the Jacobian.
482:       - Currently, all PETSc parallel matrix formats are partitioned by
483:         contiguous chunks of rows across the processors. The "grow"
484:         parameter computed below specifies the global row number
485:         corresponding to each local grid point.
486:       - Each processor needs to insert only elements that it owns
487:         locally (but any non-local elements will be sent to the
488:         appropriate processor during matrix assembly).
489:       - Always specify global row and columns of matrix entries.
490:       - Here, we set all entries for a particular row at once.
491:   */
492:   for (j=ys; j<ys+ym; j++) {
493:     row = (j - gys)*gxm + xs - gxs - 1;
494:     for (i=xs; i<xs+xm; i++) {
495:       row++;
496:       grow = ltog[row];
497:       /* boundary points */
498:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
499:         MatSetValues(jac,1,&grow,1,&grow,&one,INSERT_VALUES);
500:         continue;
501:       }
502:       /* interior grid points */
503:       v[0] = -hxdhy; col[0] = ltog[row - gxm];
504:       v[1] = -hydhx; col[1] = ltog[row - 1];
505:       v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = grow;
506:       v[3] = -hydhx; col[3] = ltog[row + 1];
507:       v[4] = -hxdhy; col[4] = ltog[row + gxm];
508:       MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
509:     }
510:   }
511:   ISLocalToGlobalMappingRestoreIndices(ltogm,&ltog);

513:   /*
514:      Assemble matrix, using the 2-step process:
515:        MatAssemblyBegin(), MatAssemblyEnd().
516:      By placing code between these two statements, computations can be
517:      done while messages are in transition.
518:   */
519:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
520:   VecRestoreArray(localX,&x);
521:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);

523:   return 0;
524: }