Actual source code: ex5s.c

petsc-3.4.5 2014-06-29
  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 mpiexec 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: */

 50: /* ------------------------------------------------------------------------

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

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

 57:     with boundary conditions

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

 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:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 74:      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: */
 94: extern int FormFunction(SNES,Vec,Vec,void*),FormInitialGuess(AppCtx*,Vec);
 95: extern int FormFunctionFortran(SNES,Vec,Vec,void*);

 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:   PetscBool      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(NULL,"-mx",&user.mx,NULL);
129:   PetscOptionsGetInt(NULL,"-my",&user.my,NULL);
130:   PetscOptionsGetReal(NULL,"-par",&user.param,NULL);
131:   if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) SETERRQ(PETSC_COMM_SELF,1,"Lambda is out of range");
132:   N = user.mx*user.my;

134:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135:      Create nonlinear solver context
136:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

138:   SNESCreate(PETSC_COMM_WORLD,&snes);

140:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141:      Create vector data structures; set function evaluation routine
142:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

153:   PetscOptionsHasName(NULL,"-use_fortran_function",&flg);
154:   if (flg) fnc = FormFunctionFortran;
155:   else     fnc = FormFunction;

157:   /*
158:      Set function evaluation routine and vector
159:   */
160:   SNESSetFunction(snes,r,fnc,&user);

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

168:     Thus this example uses the PETSc Jacobian calculations via finite differencing
169:     to approximate the Jacobian
170:   */

172:   /*

174:   */
175:   VecGetOwnershipRange(r,&rstart,&rend);
176:   PetscMalloc((rend-rstart)*sizeof(PetscInt),&colors);
177:   for (i=rstart; i<rend; i++) colors[i - rstart] = 3*((i/user.mx) % 3) + (i % 3);

179:   ISColoringCreate(PETSC_COMM_WORLD,3*2+2,rend-rstart,colors,&iscoloring);
180:   PetscFree(colors);

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

216:   /*
217:        Create the data structure that SNESComputeJacobianDefaultColor() uses
218:        to compute the actual Jacobians via finite differences.
219:   */
220:   MatFDColoringCreate(J,iscoloring,&fdcoloring);
221:   MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))fnc,&user);
222:   MatFDColoringSetFromOptions(fdcoloring);
223:   /*
224:         Tell SNES to use the routine SNESComputeJacobianDefaultColor()
225:       to compute Jacobians.
226:   */
227:   SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring);
228:   ISColoringDestroy(&iscoloring);


231:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232:      Customize nonlinear solver; set runtime options
233:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

235:   /*
236:      Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
237:   */
238:   SNESSetFromOptions(snes);

240:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241:      Evaluate initial guess; then solve nonlinear system
242:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243:   /*
244:      Note: The user should initialize the vector, x, with the initial guess
245:      for the nonlinear solver prior to calling SNESSolve().  In particular,
246:      to employ an initial guess of zero, the user should explicitly set
247:      this vector to zero by calling VecSet().
248:   */
249:   FormInitialGuess(&user,x);
250:   SNESSolve(snes,NULL,x);
251:   SNESGetIterationNumber(snes,&its);
252:   PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);

254:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
255:      Free work space.  All PETSc objects should be destroyed when they
256:      are no longer needed.
257:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
258:   VecDestroy(&x);
259:   VecDestroy(&r);
260:   SNESDestroy(&snes);
261:   PetscFinalize();

263:   return 0;
264: }
265: /* ------------------------------------------------------------------- */

269: /*
270:    FormInitialGuess - Forms initial approximation.

272:    Input Parameters:
273:    user - user-defined application context
274:    X - vector

276:    Output Parameter:
277:    X - vector
278:  */
279: int FormInitialGuess(AppCtx *user,Vec X)
280: {
281:   int         i,j,row,mx,my,ierr;
282:   PetscReal   one = 1.0,lambda,temp1,temp,hx,hy,hxdhy,hydhx,sc;
283:   PetscScalar *x;

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

300:   mx = user->mx;               my = user->my;        lambda = user->param;
301:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
302:   sc = hx*hy*lambda;        hxdhy = hx/hy;            hydhx = hy/hx;

304:   temp1 = lambda/(lambda + one);

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

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

333:   /*
334:      Restore vector
335:   */
336:   VecRestoreArray(X,&x);

338:   PetscBarrier((PetscObject)X);
339:   return 0;
340: }
341: /* ------------------------------------------------------------------- */
344: /*
345:    FormFunction - Evaluates nonlinear function, F(x).

347:    Input Parameters:
348: .  snes - the SNES context
349: .  X - input vector
350: .  ptr - optional user-defined context, as set by SNESSetFunction()

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

362:   /*
363:       Process 0 has to wait for all other processes to get here
364:    before proceeding to write in the shared vector
365:   */
366:   PetscBarrier((PetscObject)X);

368:   if (user->rank) {
369:     /*
370:        All the non-busy processors have to wait here for process 0 to finish
371:        evaluating the function; otherwise they will start using the vector values
372:        before they have been computed
373:     */
374:     PetscBarrier((PetscObject)X);
375:     return 0;
376:   }

378:   mx = user->mx;            my = user->my;            lambda = user->param;
379:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
380:   sc = hx*hy*lambda;        hxdhy = hx/hy;            hydhx = hy/hx;

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

388:   /*
389:       The next line tells the SGI compiler that x and f contain no overlapping
390:     regions and thus it can use addition optimizations.
391:   */
392: #pragma arl(4)
393: #pragma distinct (*x,*f)
394: #pragma no side effects (exp)

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

413:   /*
414:      Restore vectors
415:   */
416:   VecRestoreArray(X,&x);
417:   VecRestoreArray(F,&f);

419:   PetscLogFlops(11.0*(mx-2)*(my-2))
420:   PetscBarrier((PetscObject)X);
421:   return 0;
422: }

424: #if defined(PETSC_HAVE_FORTRAN_CAPS)
425: #define applicationfunctionfortran_ APPLICATIONFUNCTIONFORTRAN
426: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
427: #define applicationfunctionfortran_ applicationfunctionfortran
428: #endif

430: /* ------------------------------------------------------------------- */
433: /*
434:    FormFunctionFortran - Evaluates nonlinear function, F(x) in Fortran.

436: */
437: int FormFunctionFortran(SNES snes,Vec X,Vec F,void *ptr)
438: {
439:   AppCtx      *user = (AppCtx*)ptr;
440:   int         ierr;
441:   PetscScalar *x,*f;

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