Actual source code: ex1.c

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

  4: static char help[] = "Solves a nonlinear system on 1 processor with SNES. We\n\
  5: solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular domain.\n\
  6: This example also illustrates the use of matrix coloring.  Runtime options include:\n\
  7:   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
  8:      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\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\n";

 12: /*T
 13:    Concepts: SNES^sequential Bratu example
 14:    Processors: 1
 15: T*/

 17: /* ------------------------------------------------------------------------

 19:     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 20:     the partial differential equation
 21:   
 22:             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
 23:   
 24:     with boundary conditions
 25:    
 26:              u = 0  for  x = 0, x = 1, y = 0, y = 1.
 27:   
 28:     A finite difference approximation with the usual 5-point stencil
 29:     is used to discretize the boundary value problem to obtain a nonlinear 
 30:     system of equations.

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

 34:   ------------------------------------------------------------------------- */

 36: /* 
 37:    Include "petscsnes.h" so that we can use SNES solvers.  Note that
 38:    this file automatically includes:
 39:      petsc.h       - base PETSc routines   petscvec.h - vectors
 40:      petscsys.h    - system routines       petscmat.h - matrices
 41:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 42:      petscviewer.h - viewers               petscpc.h  - preconditioners
 43:      petscksp.h   - linear solvers
 44: */

 46:  #include petscsnes.h

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

 59: /* 
 60:    User-defined routines
 61: */

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

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

 85:   /*
 86:      Initialize problem parameters
 87:   */
 88:   user.mx = 4; user.my = 4; user.param = 6.0;
 89:   PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
 90:   PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
 91:   PetscOptionsGetReal(PETSC_NULL,"-par",&user.param,PETSC_NULL);
 92:   if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
 93:     SETERRQ(1,"Lambda is out of range");
 94:   }
 95:   N = user.mx*user.my;
 96: 
 97:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 98:      Create nonlinear solver context
 99:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

101:   SNESCreate(PETSC_COMM_WORLD,&snes);

103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:      Create vector data structures; set function evaluation routine
105:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

107:   VecCreate(PETSC_COMM_WORLD,&x);
108:   VecSetSizes(x,PETSC_DECIDE,N);
109:   VecSetFromOptions(x);
110:   VecDuplicate(x,&r);

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

122:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123:      Create matrix data structure; set Jacobian evaluation routine
124:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

145:   if (fd_coloring) {
146:     ISColoring   iscoloring;
147:     MatStructure str;

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

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

194:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195:      Customize nonlinear solver; set runtime options
196:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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


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

236:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237:      Free work space.  All PETSc objects should be destroyed when they
238:      are no longer needed.
239:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

241:   if (!matrix_free) {
242:     MatDestroy(J);
243:   }
244:   if (fd_coloring) {
245:     MatFDColoringDestroy(fdcoloring);
246:   }
247:   VecDestroy(x);
248:   VecDestroy(r);
249:   SNESDestroy(snes);
250:   PetscFinalize();

252:   return 0;
253: }
254: /* ------------------------------------------------------------------- */
257: /* 
258:    FormInitialGuess - Forms initial approximation.

260:    Input Parameters:
261:    user - user-defined application context
262:    X - vector

264:    Output Parameter:
265:    X - vector
266:  */
267: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
268: {
269:   PetscInt       i,j,row,mx,my;
271:   PetscReal      lambda,temp1,temp,hx,hy;
272:   PetscScalar    *x;

274:   mx         = user->mx;
275:   my         = user->my;
276:   lambda = user->param;

278:   hx    = 1.0 / (PetscReal)(mx-1);
279:   hy    = 1.0 / (PetscReal)(my-1);

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

302:   /*
303:      Restore vector
304:   */
305:   VecRestoreArray(X,&x);
306:   return 0;
307: }
308: /* ------------------------------------------------------------------- */
311: /* 
312:    FormFunction - Evaluates nonlinear function, F(x).

314:    Input Parameters:
315: .  snes - the SNES context
316: .  X - input vector
317: .  ptr - optional user-defined context, as set by SNESSetFunction()

319:    Output Parameter:
320: .  F - function vector
321:  */
322: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
323: {
324:   AppCtx         *user = (AppCtx*)ptr;
325:   PetscInt       i,j,row,mx,my;
327:   PetscReal      two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx;
328:   PetscScalar    ut,ub,ul,ur,u,uxx,uyy,sc,*x,*f;

330:   mx         = user->mx;
331:   my         = user->my;
332:   lambda = user->param;
333:   hx     = one / (PetscReal)(mx-1);
334:   hy     = one / (PetscReal)(my-1);
335:   sc     = hx*hy;
336:   hxdhy  = hx/hy;
337:   hydhx  = hy/hx;

339:   /*
340:      Get pointers to vector data
341:   */
342:   VecGetArray(X,&x);
343:   VecGetArray(F,&f);

345:   /*
346:      Compute function 
347:   */
348:   for (j=0; j<my; j++) {
349:     for (i=0; i<mx; i++) {
350:       row = i + j*mx;
351:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
352:         f[row] = x[row];
353:         continue;
354:       }
355:       u = x[row];
356:       ub = x[row - mx];
357:       ul = x[row - 1];
358:       ut = x[row + mx];
359:       ur = x[row + 1];
360:       uxx = (-ur + two*u - ul)*hydhx;
361:       uyy = (-ut + two*u - ub)*hxdhy;
362:       f[row] = uxx + uyy - sc*lambda*PetscExpScalar(u);
363:     }
364:   }

366:   /*
367:      Restore vectors
368:   */
369:   VecRestoreArray(X,&x);
370:   VecRestoreArray(F,&f);
371:   return 0;
372: }
373: /* ------------------------------------------------------------------- */
376: /*
377:    FormJacobian - Evaluates Jacobian matrix.

379:    Input Parameters:
380: .  snes - the SNES context
381: .  x - input vector
382: .  ptr - optional user-defined context, as set by SNESSetJacobian()

384:    Output Parameters:
385: .  A - Jacobian matrix
386: .  B - optionally different preconditioning matrix
387: .  flag - flag indicating matrix structure
388: */
389: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
390: {
391:   AppCtx         *user = (AppCtx*)ptr;   /* user-defined applicatin context */
392:   Mat            jac = *J;                /* Jacobian matrix */
393:   PetscInt       i,j,row,mx,my,col[5];
395:   PetscScalar    two = 2.0,one = 1.0,lambda,v[5],sc,*x;
396:   PetscReal      hx,hy,hxdhy,hydhx;

398:   mx         = user->mx;
399:   my         = user->my;
400:   lambda = user->param;
401:   hx     = 1.0 / (PetscReal)(mx-1);
402:   hy     = 1.0 / (PetscReal)(my-1);
403:   sc     = hx*hy;
404:   hxdhy  = hx/hy;
405:   hydhx  = hy/hx;

407:   /*
408:      Get pointer to vector data
409:   */
410:   VecGetArray(X,&x);

412:   /* 
413:      Compute entries of the Jacobian
414:   */
415:   for (j=0; j<my; j++) {
416:     for (i=0; i<mx; i++) {
417:       row = i + j*mx;
418:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
419:         MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES);
420:         continue;
421:       }
422:       v[0] = -hxdhy; col[0] = row - mx;
423:       v[1] = -hydhx; col[1] = row - 1;
424:       v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = row;
425:       v[3] = -hydhx; col[3] = row + 1;
426:       v[4] = -hxdhy; col[4] = row + mx;
427:       MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES);
428:     }
429:   }

431:   /*
432:      Restore vector
433:   */
434:   VecRestoreArray(X,&x);

436:   /* 
437:      Assemble matrix
438:   */
439:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
440:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);

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