Actual source code: ex1.c

  1: /*$Id: ex1.c,v 1.87 2001/08/27 17:31:28 bsmith Exp $*/

  3: /* Program usage:  ex4 [-help] [all PETSc options] */

  5: static char help[] = "Solves a nonlinear system on 1 processor with SNES. Wen
  6: solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular domain.n
  7: This example also illustrates the use of matrix coloring.  Runtime options include:n
  8:   -par <parameter>, where <parameter> indicates the problem's nonlinearityn
  9:      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)n
 10:   -mx <xg>, where <xg> = number of grid points in the x-directionn
 11:   -my <yg>, where <yg> = number of grid points in the y-directionnn";

 13: /*T
 14:    Concepts: SNES^sequential Bratu example
 15:    Processors: 1
 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.
 28:   
 29:     A finite difference approximation with the usual 5-point stencil
 30:     is used to discretize the boundary value problem to obtain a nonlinear 
 31:     system of equations.

 33:     The parallel version of this code is snes/examples/tutorials/ex5.c

 35:   ------------------------------------------------------------------------- */

 37: /* 
 38:    Include "petscsnes.h" so that we can use SNES solvers.  Note that
 39:    this 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: */

 47:  #include petscsnes.h

 49: /* 
 50:    User-defined application context - contains data needed by the 
 51:    application-provided call-back routines, FormJacobian() and
 52:    FormFunction().
 53: */
 54: typedef struct {
 55:       PetscReal   param;        /* test problem parameter */
 56:       int         mx;           /* Discretization in x-direction */
 57:       int         my;           /* Discretization in y-direction */
 58: } AppCtx;

 60: /* 
 61:    User-defined routines
 62: */
 63: extern int FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
 64: extern int FormFunction(SNES,Vec,Vec,void*);
 65: extern int FormInitialGuess(AppCtx*,Vec);

 67: int main(int argc,char **argv)
 68: {
 69:   SNES           snes;                 /* nonlinear solver context */
 70:   Vec            x,r;                 /* solution, residual vectors */
 71:   Mat            J;                    /* Jacobian matrix */
 72:   AppCtx         user;                 /* user-defined application context */
 73:   int            i,ierr,its,N,size,hist_its[50];
 74:   PetscReal      bratu_lambda_max = 6.81,bratu_lambda_min = 0.,history[50];
 75:   MatFDColoring  fdcoloring;
 76:   PetscTruth     matrix_free,flg,fd_coloring;

 78:   PetscInitialize(&argc,&argv,(char *)0,help);
 79:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 80:   if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");

 82:   /*
 83:      Initialize problem parameters
 84:   */
 85:   user.mx = 4; user.my = 4; user.param = 6.0;
 86:   PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
 87:   PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
 88:   PetscOptionsGetReal(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:   }
 92:   N = user.mx*user.my;
 93: 
 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:      Create nonlinear solver context
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 98:   SNESCreate(PETSC_COMM_WORLD,&snes);

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:      Create vector data structures; set function evaluation routine
102:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

104:   VecCreate(PETSC_COMM_WORLD,&x);
105:   VecSetSizes(x,PETSC_DECIDE,N);
106:   VecSetFromOptions(x);
107:   VecDuplicate(x,&r);

109:   /* 
110:      Set function evaluation routine and vector.  Whenever the nonlinear
111:      solver needs to evaluate the nonlinear function, it will call this
112:      routine.
113:       - Note that the final routine argument is the user-defined
114:         context that provides application-specific data for the
115:         function evaluation routine.
116:   */
117:   SNESSetFunction(snes,r,FormFunction,(void *)&user);

119:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120:      Create matrix data structure; set Jacobian evaluation routine
121:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

123:   /*
124:      Create matrix. Here we only approximately preallocate storage space
125:      for the Jacobian.  See the users manual for a discussion of better 
126:      techniques for preallocating matrix memory.
127:   */
128:   PetscOptionsHasName(PETSC_NULL,"-snes_mf",&matrix_free);
129:   if (!matrix_free) {
130:     MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,PETSC_NULL,&J);
131:   }

133:   /*
134:      This option will cause the Jacobian to be computed via finite differences
135:     efficiently using a coloring of the columns of the matrix.
136:   */
137:   PetscOptionsHasName(PETSC_NULL,"-snes_fd_coloring",&fd_coloring);

139:   if (matrix_free && fd_coloring)  SETERRQ(1,"Use only one of -snes_mf, -snes_fd_coloring options!n
140:                                                 You can do -snes_mf_operator -snes_fd_coloring");

142:   if (fd_coloring) {
143:     ISColoring   iscoloring;
144:     MatStructure str;

146:     /* 
147:       This initializes the nonzero structure of the Jacobian. This is artificial
148:       because clearly if we had a routine to compute the Jacobian we won't need
149:       to use finite differences.
150:     */
151:     FormJacobian(snes,x,&J,&J,&str,&user);

153:     /*
154:        Color the matrix, i.e. determine groups of columns that share no common 
155:       rows. These columns in the Jacobian can all be computed simulataneously.
156:     */
157:     MatGetColoring(J,MATCOLORING_NATURAL,&iscoloring);
158:     /*
159:        Create the data structure that SNESDefaultComputeJacobianColor() uses
160:        to compute the actual Jacobians via finite differences.
161:     */
162:     MatFDColoringCreate(J,iscoloring,&fdcoloring);
163:     MatFDColoringSetFunction(fdcoloring,(int (*)(void))FormFunction,&user);
164:     MatFDColoringSetFromOptions(fdcoloring);
165:     /*
166:         Tell SNES to use the routine SNESDefaultComputeJacobianColor()
167:       to compute Jacobians.
168:     */
169:     SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,fdcoloring);
170:     ISColoringDestroy(iscoloring);
171:   }
172:   /* 
173:      Set Jacobian matrix data structure and default Jacobian evaluation
174:      routine.  Whenever the nonlinear solver needs to compute the
175:      Jacobian matrix, it will call this routine.
176:       - Note that the final routine argument is the user-defined
177:         context that provides application-specific data for the
178:         Jacobian evaluation routine.
179:       - The user can override with:
180:          -snes_fd : default finite differencing approximation of Jacobian
181:          -snes_mf : matrix-free Newton-Krylov method with no preconditioning
182:                     (unless user explicitly sets preconditioner) 
183:          -snes_mf_operator : form preconditioning matrix as set by the user,
184:                              but use matrix-free approx for Jacobian-vector
185:                              products within Newton-Krylov method
186:   */
187:   else if (!matrix_free) {
188:     SNESSetJacobian(snes,J,J,FormJacobian,(void*)&user);
189:   }

191:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192:      Customize nonlinear solver; set runtime options
193:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

200:   /*
201:      Set array that saves the function norms.  This array is intended
202:      when the user wants to save the convergence history for later use
203:      rather than just to view the function norms via -snes_monitor.
204:   */
205:   SNESSetConvergenceHistory(snes,history,hist_its,50,PETSC_TRUE);

207:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208:      Evaluate initial guess; then solve nonlinear system
209:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
210:   /*
211:      Note: The user should initialize the vector, x, with the initial guess
212:      for the nonlinear solver prior to calling SNESSolve().  In particular,
213:      to employ an initial guess of zero, the user should explicitly set
214:      this vector to zero by calling VecSet().
215:   */
216:   FormInitialGuess(&user,x);
217:   SNESSolve(snes,x,&its);
218:   PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %dn",its);


221:   /* 
222:      Print the convergence history.  This is intended just to demonstrate
223:      use of the data attained via SNESSetConvergenceHistory().  
224:   */
225:   PetscOptionsHasName(PETSC_NULL,"-print_history",&flg);
226:   if (flg) {
227:     for (i=0; i<its+1; i++) {
228:       PetscPrintf(PETSC_COMM_WORLD,"iteration %d: Linear iterations %d Function norm = %gn",i,hist_its[i],history[i]);
229:     }
230:   }

232:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233:      Free work space.  All PETSc objects should be destroyed when they
234:      are no longer needed.
235:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

237:   if (!matrix_free) {
238:     MatDestroy(J);
239:   }
240:   if (fd_coloring) {
241:     MatFDColoringDestroy(fdcoloring);
242:   }
243:   VecDestroy(x);
244:   VecDestroy(r);
245:   SNESDestroy(snes);
246:   PetscFinalize();

248:   return 0;
249: }
250: /* ------------------------------------------------------------------- */
251: /* 
252:    FormInitialGuess - Forms initial approximation.

254:    Input Parameters:
255:    user - user-defined application context
256:    X - vector

258:    Output Parameter:
259:    X - vector
260:  */
261: int FormInitialGuess(AppCtx *user,Vec X)
262: {
263:   int         i,j,row,mx,my,ierr;
264:   PetscReal   lambda,temp1,temp,hx,hy;
265:   PetscScalar *x;

267:   mx         = user->mx;
268:   my         = user->my;
269:   lambda = user->param;

271:   hx    = 1.0 / (PetscReal)(mx-1);
272:   hy    = 1.0 / (PetscReal)(my-1);

274:   /*
275:      Get a pointer to vector data.
276:        - For default PETSc vectors, VecGetArray() returns a pointer to
277:          the data array.  Otherwise, the routine is implementation dependent.
278:        - You MUST call VecRestoreArray() when you no longer need access to
279:          the array.
280:   */
281:   VecGetArray(X,&x);
282:   temp1 = lambda/(lambda + 1.0);
283:   for (j=0; j<my; j++) {
284:     temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
285:     for (i=0; i<mx; i++) {
286:       row = i + j*mx;
287:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
288:         x[row] = 0.0;
289:         continue;
290:       }
291:       x[row] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
292:     }
293:   }

295:   /*
296:      Restore vector
297:   */
298:   VecRestoreArray(X,&x);
299:   return 0;
300: }
301: /* ------------------------------------------------------------------- */
302: /* 
303:    FormFunction - Evaluates nonlinear function, F(x).

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

310:    Output Parameter:
311: .  F - function vector
312:  */
313: int FormFunction(SNES snes,Vec X,Vec F,void *ptr)
314: {
315:   AppCtx       *user = (AppCtx*)ptr;
316:   int          ierr,i,j,row,mx,my;
317:   PetscReal    two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx;
318:   PetscScalar  ut,ub,ul,ur,u,uxx,uyy,sc,*x,*f;

320:   mx         = user->mx;
321:   my         = user->my;
322:   lambda = user->param;
323:   hx     = one / (PetscReal)(mx-1);
324:   hy     = one / (PetscReal)(my-1);
325:   sc     = hx*hy;
326:   hxdhy  = hx/hy;
327:   hydhx  = hy/hx;

329:   /*
330:      Get pointers to vector data
331:   */
332:   VecGetArray(X,&x);
333:   VecGetArray(F,&f);

335:   /*
336:      Compute function 
337:   */
338:   for (j=0; j<my; j++) {
339:     for (i=0; i<mx; i++) {
340:       row = i + j*mx;
341:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
342:         f[row] = x[row];
343:         continue;
344:       }
345:       u = x[row];
346:       ub = x[row - mx];
347:       ul = x[row - 1];
348:       ut = x[row + mx];
349:       ur = x[row + 1];
350:       uxx = (-ur + two*u - ul)*hydhx;
351:       uyy = (-ut + two*u - ub)*hxdhy;
352:       f[row] = uxx + uyy - sc*lambda*PetscExpScalar(u);
353:     }
354:   }

356:   /*
357:      Restore vectors
358:   */
359:   VecRestoreArray(X,&x);
360:   VecRestoreArray(F,&f);
361:   return 0;
362: }
363: /* ------------------------------------------------------------------- */
364: /*
365:    FormJacobian - Evaluates Jacobian matrix.

367:    Input Parameters:
368: .  snes - the SNES context
369: .  x - input vector
370: .  ptr - optional user-defined context, as set by SNESSetJacobian()

372:    Output Parameters:
373: .  A - Jacobian matrix
374: .  B - optionally different preconditioning matrix
375: .  flag - flag indicating matrix structure
376: */
377: int FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
378: {
379:   AppCtx      *user = (AppCtx*)ptr;   /* user-defined applicatin context */
380:   Mat         jac = *J;                /* Jacobian matrix */
381:   int         i,j,row,mx,my,col[5],ierr;
382:   PetscScalar two = 2.0,one = 1.0,lambda,v[5],sc,*x;
383:   PetscReal   hx,hy,hxdhy,hydhx;

385:   mx         = user->mx;
386:   my         = user->my;
387:   lambda = user->param;
388:   hx     = 1.0 / (PetscReal)(mx-1);
389:   hy     = 1.0 / (PetscReal)(my-1);
390:   sc     = hx*hy;
391:   hxdhy  = hx/hy;
392:   hydhx  = hy/hx;

394:   /*
395:      Get pointer to vector data
396:   */
397:   VecGetArray(X,&x);

399:   /* 
400:      Compute entries of the Jacobian
401:   */
402:   for (j=0; j<my; j++) {
403:     for (i=0; i<mx; i++) {
404:       row = i + j*mx;
405:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
406:         MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES);
407:         continue;
408:       }
409:       v[0] = -hxdhy; col[0] = row - mx;
410:       v[1] = -hydhx; col[1] = row - 1;
411:       v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = row;
412:       v[3] = -hydhx; col[3] = row + 1;
413:       v[4] = -hxdhy; col[4] = row + mx;
414:       MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES);
415:     }
416:   }

418:   /*
419:      Restore vector
420:   */
421:   VecRestoreArray(X,&x);

423:   /* 
424:      Assemble matrix
425:   */
426:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
427:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);

429:   /*
430:      Set flag to indicate that the Jacobian matrix retains an identical
431:      nonzero structure throughout all nonlinear iterations (although the
432:      values of the entries change). Thus, we can save some work in setting
433:      up the preconditioner (e.g., no need to redo symbolic factorization for
434:      ILU/ICC preconditioners).
435:       - If the nonzero structure of the matrix is different during
436:         successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
437:         must be used instead.  If you are unsure whether the matrix
438:         structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
439:       - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
440:         believes your assertion and does not check the structure
441:         of the matrix.  If you erroneously claim that the structure
442:         is the same when it actually is not, the new preconditioner
443:         will not function correctly.  Thus, use this optimization
444:         feature with caution!
445:   */
446:   *flag = SAME_NONZERO_PATTERN;
447:   return 0;
448: }