Actual source code: ex5s.c

  2: static char help[] = "2d Bratu problem in shared memory parallel with SNES.\n\
  3: We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\
  4: domain, uses SHARED MEMORY to evaluate the user function.\n\
  5: The command line options include:\n\
  6:   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
  7:      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\
  8:   -mx <xg>, where <xg> = number of grid points in the x-direction\n\
  9:   -my <yg>, where <yg> = number of grid points in the y-direction\n\
 10:   -use_fortran_function: use Fortran coded function, rather than C\n";

 12: /*
 13:              This code compiles ONLY on SGI systems
 14:             ========================================
 15: */
 16: /*T
 17:    Concepts: SNES^parallel Bratu example
 18:    Concepts: shared memory
 19:    Processors: n
 20: T*/

 22: /*

 24:      Programming model: Combination of 
 25:         1) MPI message passing for PETSc routines
 26:         2) automatic loop parallism (using shared memory) for user
 27:            provided function.

 29:        While the user function is being evaluated all MPI processes except process
 30:      0 blocks. Process zero spawns nt threads to evaluate the user function. Once 
 31:      the user function is complete, the worker threads are suspended and all the MPI processes
 32:      continue.

 34:      Other useful options:

 36:        -snes_mf : use matrix free operator and no preconditioner
 37:        -snes_mf_operator : use matrix free operator but compute Jacobian via
 38:                            finite differences to form preconditioner

 40:        Environmental variable:

 42:          setenv MPC_NUM_THREADS nt <- set number of threads processor 0 should 
 43:                                       use to evaluate user provided function

 45:        Note: The number of MPI processes (set with the mpirun option -np) can 
 46:        be set completely independently from the number of threads process 0 
 47:        uses to evaluate the function (though usually one would make them the same).
 48: */
 49: 
 50: /* ------------------------------------------------------------------------

 52:     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 53:     the partial differential equation
 54:   
 55:             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
 56:   
 57:     with boundary conditions
 58:    
 59:              u = 0  for  x = 0, x = 1, y = 0, y = 1.
 60:   
 61:     A finite difference approximation with the usual 5-point stencil
 62:     is used to discretize the boundary value problem to obtain a nonlinear 
 63:     system of equations.

 65:     The uniprocessor version of this code is snes/examples/tutorials/ex4.c
 66:     A parallel distributed memory version is snes/examples/tutorials/ex5.c and ex5f.F

 68:   ------------------------------------------------------------------------- */

 70: /* 
 71:    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
 72:    file automatically includes:
 73:      petsc.h       - base PETSc routines   petscvec.h - vectors
 74:      petscsys.h    - system routines       petscmat.h - matrices
 75:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 76:      petscviewer.h - viewers               petscpc.h  - preconditioners
 77:      petscksp.h   - linear solvers
 78: */
 79:  #include petscsnes.h

 81: /* 
 82:    User-defined application context - contains data needed by the 
 83:    application-provided call-back routines   FormFunction().
 84: */
 85: typedef struct {
 86:    PetscReal   param;          /* test problem parameter */
 87:    int         mx,my;          /* discretization in x, y directions */
 88:    int         rank;           /* processor rank */
 89: } AppCtx;

 91: /* 
 92:    User-defined routines
 93: */

 99: /* 
100:     The main program is written in C while the user provided function
101:  is given in both Fortran and C. The main program could also be written 
102:  in Fortran; the ONE PROBLEM is that VecGetArray() cannot be called from 
103:  Fortran on the SGI machines; thus the routine FormFunctionFortran() must
104:  be written in C.
105: */
106: int main(int argc,char **argv)
107: {
108:   SNES           snes;                /* nonlinear solver */
109:   Vec            x,r;                 /* solution, residual vectors */
110:   AppCtx         user;                /* user-defined work context */
111:   int            its;                 /* iterations for convergence */
112:   int            N,ierr,rstart,rend,*colors,i,ii,ri,rj;
113:   PetscErrorCode (*fnc)(SNES,Vec,Vec,void*);
114:   PetscReal      bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
115:   MatFDColoring  fdcoloring;
116:   ISColoring     iscoloring;
117:   Mat            J;
118:   PetscScalar    zero = 0.0;
119:   PetscTruth     flg;

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

124:   /*
125:      Initialize problem parameters
126:   */
127:   user.mx = 4; user.my = 4; user.param = 6.0;
128:   PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
129:   PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
130:   PetscOptionsGetReal(PETSC_NULL,"-par",&user.param,PETSC_NULL);
131:   if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
132:     SETERRQ(1,"Lambda is out of range");
133:   }
134:   N = user.mx*user.my;

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:      Create nonlinear solver context
138:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

140:   SNESCreate(PETSC_COMM_WORLD,&snes);

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:      Create vector data structures; set function evaluation routine
144:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

146:   /*
147:       The routine VecCreateShared() creates a parallel vector with each processor
148:     assigned its own segment, BUT, in addition, the first processor has access to the 
149:     entire array. This is to allow the users function to be based on loop level
150:     parallelism rather than MPI.
151:   */
152:   VecCreateShared(PETSC_COMM_WORLD,PETSC_DECIDE,N,&x);
153:   VecDuplicate(x,&r);

155:   PetscOptionsHasName(PETSC_NULL,"-use_fortran_function",&flg);
156:   if (flg) {
157:     fnc = FormFunctionFortran;
158:   } else {
159:     fnc = FormFunction;
160:   }

162:   /* 
163:      Set function evaluation routine and vector
164:   */
165:   SNESSetFunction(snes,r,fnc,&user);

167:   /*
168:        Currently when using VecCreateShared() and using loop level parallelism
169:     to automatically parallelise the user function it makes no sense for the 
170:     Jacobian to be computed via loop level parallelism, because all the threads
171:     would be simultaneously calling MatSetValues() causing a bottle-neck.

173:     Thus this example uses the PETSc Jacobian calculations via finite differencing
174:     to approximate the Jacobian
175:   */

177:   /*

179:   */
180:   VecGetOwnershipRange(r,&rstart,&rend);
181:   PetscMalloc((rend-rstart)*sizeof(PetscInt),&colors);
182:   for (i=rstart; i<rend; i++) {
183:     colors[i - rstart] = 3*((i/user.mx) % 3) + (i % 3);
184:   }
185:   ISColoringCreate(PETSC_COMM_WORLD,3*2+2,rend-rstart,colors,&iscoloring);
186:   PetscFree(colors);

188:   /*
189:      Create and set the nonzero pattern for the Jacobian: This is not done 
190:      particularly efficiently. One should process the boundary nodes separately and 
191:      then use a simple loop for the interior nodes.
192:        Note that for this code we use the "natural" number of the nodes on the 
193:      grid (since that is what is good for the user provided function). In the 
194:      DA examples we must use the DA numbering where each processor is assigned a
195:      chunk of data.
196:   */
197:   MatCreateMPIAIJ(PETSC_COMM_WORLD,rend-rstart,rend-rstart,N,
198:                          N,5,0,0,0,&J);
199:   for (i=rstart; i<rend; i++) {
200:     rj = i % user.mx;         /* column in grid */
201:     ri = i / user.mx;         /* row in grid */
202:     if (ri != 0) {     /* first row does not have neighbor below */
203:       ii   = i - user.mx;
204:       MatSetValues(J,1,&i,1,&ii,&zero,INSERT_VALUES);
205:     }
206:     if (ri != user.my - 1) { /* last row does not have neighbors above */
207:       ii   = i + user.mx;
208:       MatSetValues(J,1,&i,1,&ii,&zero,INSERT_VALUES);
209:     }
210:     if (rj != 0) {     /* first column does not have neighbor to left */
211:       ii   = i - 1;
212:       MatSetValues(J,1,&i,1,&ii,&zero,INSERT_VALUES);
213:     }
214:     if (rj != user.mx - 1) {     /* last column does not have neighbor to right */
215:       ii   = i + 1;
216:       MatSetValues(J,1,&i,1,&ii,&zero,INSERT_VALUES);
217:     }
218:     MatSetValues(J,1,&i,1,&i,&zero,INSERT_VALUES);
219:   }
220:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
221:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);

223:   /*
224:        Create the data structure that SNESDefaultComputeJacobianColor() uses
225:        to compute the actual Jacobians via finite differences.
226:   */
227:   MatFDColoringCreate(J,iscoloring,&fdcoloring);
228:   MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))fnc,&user);
229:   MatFDColoringSetFromOptions(fdcoloring);
230:   /*
231:         Tell SNES to use the routine SNESDefaultComputeJacobianColor()
232:       to compute Jacobians.
233:   */
234:   SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,fdcoloring);
235:   ISColoringDestroy(iscoloring);


238:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239:      Customize nonlinear solver; set runtime options
240:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

242:   /*
243:      Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
244:   */
245:   SNESSetFromOptions(snes);

247:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248:      Evaluate initial guess; then solve nonlinear system
249:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
250:   /*
251:      Note: The user should initialize the vector, x, with the initial guess
252:      for the nonlinear solver prior to calling SNESSolve().  In particular,
253:      to employ an initial guess of zero, the user should explicitly set
254:      this vector to zero by calling VecSet().
255:   */
256:   FormInitialGuess(&user,x);
257:   SNESSolve(snes,PETSC_NULL,x);
258:   SNESGetIterationNumber(snes,&its);
259:   PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %D\n",its);

261:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262:      Free work space.  All PETSc objects should be destroyed when they
263:      are no longer needed.
264:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
265:   VecDestroy(x);
266:   VecDestroy(r);
267:   SNESDestroy(snes);
268:   PetscFinalize();

270:   return 0;
271: }
272: /* ------------------------------------------------------------------- */

276: /* 
277:    FormInitialGuess - Forms initial approximation.

279:    Input Parameters:
280:    user - user-defined application context
281:    X - vector

283:    Output Parameter:
284:    X - vector
285:  */
286: int FormInitialGuess(AppCtx *user,Vec X)
287: {
288:   int          i,j,row,mx,my,ierr;
289:   PetscReal    one = 1.0,lambda,temp1,temp,hx,hy,hxdhy,hydhx,sc;
290:   PetscScalar  *x;

292:   /*
293:       Process 0 has to wait for all other processes to get here 
294:    before proceeding to write in the shared vector
295:   */
296:   PetscBarrier((PetscObject)X);
297:   if (user->rank) {
298:      /*
299:         All the non-busy processors have to wait here for process 0 to finish
300:         evaluating the function; otherwise they will start using the vector values
301:         before they have been computed
302:      */
303:      PetscBarrier((PetscObject)X);
304:      return 0;
305:   }

307:   mx = user->mx;            my = user->my;            lambda = user->param;
308:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
309:   sc = hx*hy*lambda;        hxdhy = hx/hy;            hydhx = hy/hx;
310:   temp1 = lambda/(lambda + one);

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

321:   /*
322:      Compute initial guess over the locally owned part of the grid
323:   */
324: #pragma arl(4)
325: #pragma distinct (*x,*f)
326: #pragma no side effects (sqrt)
327:   for (j=0; j<my; j++) {
328:     temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
329:     for (i=0; i<mx; i++) {
330:       row = i + j*mx;
331:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
332:         x[row] = 0.0;
333:         continue;
334:       }
335:       x[row] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
336:     }
337:   }

339:   /*
340:      Restore vector
341:   */
342:   VecRestoreArray(X,&x);

344:   PetscBarrier((PetscObject)X);
345:   return 0;
346: }
347: /* ------------------------------------------------------------------- */
350: /* 
351:    FormFunction - Evaluates nonlinear function, F(x).

353:    Input Parameters:
354: .  snes - the SNES context
355: .  X - input vector
356: .  ptr - optional user-defined context, as set by SNESSetFunction()

358:    Output Parameter:
359: .  F - function vector
360:  */
361: int FormFunction(SNES snes,Vec X,Vec F,void *ptr)
362: {
363:   AppCtx       *user = (AppCtx*)ptr;
364:   int          ierr,i,j,row,mx,my;
365:   PetscReal    two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx,sc;
366:   PetscScalar  u,uxx,uyy,*x,*f;

368:   /*
369:       Process 0 has to wait for all other processes to get here 
370:    before proceeding to write in the shared vector
371:   */
372:   PetscBarrier((PetscObject)X);

374:   if (user->rank) {
375:      /*
376:         All the non-busy processors have to wait here for process 0 to finish
377:         evaluating the function; otherwise they will start using the vector values
378:         before they have been computed
379:      */
380:      PetscBarrier((PetscObject)X);
381:      return 0;
382:   }

384:   mx = user->mx;            my = user->my;            lambda = user->param;
385:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
386:   sc = hx*hy*lambda;        hxdhy = hx/hy;            hydhx = hy/hx;

388:   /*
389:      Get pointers to vector data
390:   */
391:   VecGetArray(X,&x);
392:   VecGetArray(F,&f);

394:   /*
395:       The next line tells the SGI compiler that x and f contain no overlapping 
396:     regions and thus it can use addition optimizations.
397:   */
398: #pragma arl(4)
399: #pragma distinct (*x,*f)
400: #pragma no side effects (exp)

402:   /*
403:      Compute function over the entire  grid
404:   */
405:   for (j=0; j<my; j++) {
406:     for (i=0; i<mx; i++) {
407:       row = i + j*mx;
408:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
409:         f[row] = x[row];
410:         continue;
411:       }
412:       u = x[row];
413:       uxx = (two*u - x[row-1] - x[row+1])*hydhx;
414:       uyy = (two*u - x[row-mx] - x[row+mx])*hxdhy;
415:       f[row] = uxx + uyy - sc*exp(u);
416:     }
417:   }

419:   /*
420:      Restore vectors
421:   */
422:   VecRestoreArray(X,&x);
423:   VecRestoreArray(F,&f);

425:   PetscLogFlops(11*(mx-2)*(my-2))
426:   PetscBarrier((PetscObject)X);
427:   return 0;
428: }

430: #if defined(PETSC_HAVE_FORTRAN_CAPS)
431: #define applicationfunctionfortran_ APPLICATIONFUNCTIONFORTRAN
432: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
433: #define applicationfunctionfortran_ applicationfunctionfortran
434: #endif

436: /* ------------------------------------------------------------------- */
439: /* 
440:    FormFunctionFortran - Evaluates nonlinear function, F(x) in Fortran.

442: */
443: int FormFunctionFortran(SNES snes,Vec X,Vec F,void *ptr)
444: {
445:   AppCtx  *user = (AppCtx*)ptr;
446:   int     ierr;
447:   PetscScalar  *x,*f;

449:   /*
450:       Process 0 has to wait for all other processes to get here 
451:    before proceeding to write in the shared vector
452:   */
453:   PetscBarrier((PetscObject)snes);
454:   if (!user->rank) {
455:     VecGetArray(X,&x);
456:     VecGetArray(F,&f);
457:     applicationfunctionfortran_(&user->param,&user->mx,&user->my,x,f,&ierr);
458:     VecRestoreArray(X,&x);
459:     VecRestoreArray(F,&f);
460:     PetscLogFlops(11*(user->mx-2)*(user->my-2))
461:   }
462:   /*
463:       All the non-busy processors have to wait here for process 0 to finish
464:       evaluating the function; otherwise they will start using the vector values
465:       before they have been computed
466:   */
467:   PetscBarrier((PetscObject)snes);
468:   return 0;
469: }