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,NULL,"-mx",&user.mx,NULL);
129: PetscOptionsGetInt(NULL,NULL,"-my",&user.my,NULL);
130: PetscOptionsGetReal(NULL,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,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: PetscMalloc1(rend-rstart,&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: MatFDColoringSetUp(J,iscoloring,fdcoloring);
224: /*
225: Tell SNES to use the routine SNESComputeJacobianDefaultColor()
226: to compute Jacobians.
227: */
228: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring);
229: ISColoringDestroy(&iscoloring);
232: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233: Customize nonlinear solver; set runtime options
234: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236: /*
237: Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
238: */
239: SNESSetFromOptions(snes);
241: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
242: Evaluate initial guess; then solve nonlinear system
243: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
244: /*
245: Note: The user should initialize the vector, x, with the initial guess
246: for the nonlinear solver prior to calling SNESSolve(). In particular,
247: to employ an initial guess of zero, the user should explicitly set
248: this vector to zero by calling VecSet().
249: */
250: FormInitialGuess(&user,x);
251: SNESSolve(snes,NULL,x);
252: SNESGetIterationNumber(snes,&its);
253: PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);
255: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
256: Free work space. All PETSc objects should be destroyed when they
257: are no longer needed.
258: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
259: VecDestroy(&x);
260: VecDestroy(&r);
261: SNESDestroy(&snes);
262: PetscFinalize();
264: return 0;
265: }
266: /* ------------------------------------------------------------------- */
270: /*
271: FormInitialGuess - Forms initial approximation.
273: Input Parameters:
274: user - user-defined application context
275: X - vector
277: Output Parameter:
278: X - vector
279: */
280: int FormInitialGuess(AppCtx *user,Vec X)281: {
282: int i,j,row,mx,my,ierr;
283: PetscReal one = 1.0,lambda,temp1,temp,hx,hy,hxdhy,hydhx,sc;
284: PetscScalar *x;
286: /*
287: Process 0 has to wait for all other processes to get here
288: before proceeding to write in the shared vector
289: */
290: PetscBarrier((PetscObject)X);
291: if (user->rank) {
292: /*
293: All the non-busy processors have to wait here for process 0 to finish
294: evaluating the function; otherwise they will start using the vector values
295: before they have been computed
296: */
297: PetscBarrier((PetscObject)X);
298: return 0;
299: }
301: mx = user->mx; my = user->my; lambda = user->param;
302: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
303: sc = hx*hy*lambda; hxdhy = hx/hy; hydhx = hy/hx;
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: Compute initial guess over the locally owned part of the grid
318: */
319: #pragma arl(4)
320: #pragma distinct (*x,*f)
321: #pragma no side effects (sqrt)
322: for (j=0; j<my; j++) {
323: temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
324: for (i=0; i<mx; i++) {
325: row = i + j*mx;
326: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
327: x[row] = 0.0;
328: continue;
329: }
330: x[row] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
331: }
332: }
334: /*
335: Restore vector
336: */
337: VecRestoreArray(X,&x);
339: PetscBarrier((PetscObject)X);
340: return 0;
341: }
342: /* ------------------------------------------------------------------- */
345: /*
346: FormFunction - Evaluates nonlinear function, F(x).
348: Input Parameters:
349: . snes - the SNES context
350: . X - input vector
351: . ptr - optional user-defined context, as set by SNESSetFunction()
353: Output Parameter:
354: . F - function vector
355: */
356: int FormFunction(SNES snes,Vec X,Vec F,void *ptr)357: {
358: AppCtx *user = (AppCtx*)ptr;
359: int ierr,i,j,row,mx,my;
360: PetscReal two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx,sc;
361: PetscScalar u,uxx,uyy,*f;
362: const PetscScalar *x;
364: /*
365: Process 0 has to wait for all other processes to get here
366: before proceeding to write in the shared vector
367: */
368: PetscBarrier((PetscObject)X);
370: if (user->rank) {
371: /*
372: All the non-busy processors have to wait here for process 0 to finish
373: evaluating the function; otherwise they will start using the vector values
374: before they have been computed
375: */
376: PetscBarrier((PetscObject)X);
377: return 0;
378: }
380: mx = user->mx; my = user->my; lambda = user->param;
381: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
382: sc = hx*hy*lambda; hxdhy = hx/hy; hydhx = hy/hx;
384: /*
385: Get pointers to vector data
386: */
387: VecGetArrayRead(X,&x);
388: VecGetArray(F,&f);
390: /*
391: The next line tells the SGI compiler that x and f contain no overlapping
392: regions and thus it can use addition optimizations.
393: */
394: #pragma arl(4)
395: #pragma distinct (*x,*f)
396: #pragma no side effects (exp)
398: /*
399: Compute function over the entire grid
400: */
401: for (j=0; j<my; j++) {
402: for (i=0; i<mx; i++) {
403: row = i + j*mx;
404: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
405: f[row] = x[row];
406: continue;
407: }
408: u = x[row];
409: uxx = (two*u - x[row-1] - x[row+1])*hydhx;
410: uyy = (two*u - x[row-mx] - x[row+mx])*hxdhy;
411: f[row] = uxx + uyy - sc*PetscExpScalar(u);
412: }
413: }
415: /*
416: Restore vectors
417: */
418: VecRestoreArrayRead(X,&x);
419: VecRestoreArray(F,&f);
421: PetscLogFlops(11.0*(mx-2)*(my-2))
422: PetscBarrier((PetscObject)X);
423: return 0;
424: }
426: #if defined(PETSC_HAVE_FORTRAN_CAPS)
427: #define applicationfunctionfortran_ APPLICATIONFUNCTIONFORTRAN428: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
429: #define applicationfunctionfortran_ applicationfunctionfortran430: #endif
432: /* ------------------------------------------------------------------- */
435: /*
436: FormFunctionFortran - Evaluates nonlinear function, F(x) in Fortran.
438: */
439: int FormFunctionFortran(SNES snes,Vec X,Vec F,void *ptr)440: {
441: AppCtx *user = (AppCtx*)ptr;
442: int ierr;
443: PetscScalar *f;
444: const PetscScalar *x;
446: /*
447: Process 0 has to wait for all other processes to get here
448: before proceeding to write in the shared vector
449: */
450: PetscBarrier((PetscObject)snes);
451: if (!user->rank) {
452: VecGetArrayRead(X,&x);
453: VecGetArray(F,&f);
454: applicationfunctionfortran_(&user->param,&user->mx,&user->my,x,f,&ierr);
455: VecRestoreArrayRead(X,&x);
456: VecRestoreArray(F,&f);
457: PetscLogFlops(11.0*(user->mx-2)*(user->my-2))
458: }
459: /*
460: All the non-busy processors have to wait here for process 0 to finish
461: evaluating the function; otherwise they will start using the vector values
462: before they have been computed
463: */
464: PetscBarrier((PetscObject)snes);
465: return 0;
466: }