Actual source code: ex5.c
1: /*$Id: ex5.c,v 1.29 2001/08/07 21:31:12 bsmith Exp $*/
3: static char help[] = "Solves a nonlinear system in parallel with SNES.\n\
4: We solve the modified Bratu problem in a 2D rectangular domain,\n\
5: using distributed arrays (DAs) to partition the parallel grid.\n\
6: The command line options include:\n\
7: -lambda <parameter>, where <parameter> indicates the problem's nonlinearity\n\
8: -kappa <parameter>, where <parameter> indicates the problem's nonlinearity\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\
11: -Nx <npx>, where <npx> = number of processors in the x-direction\n\
12: -Ny <npy>, where <npy> = number of processors in the y-direction\n\n";
14: /*T
15: Concepts: SNES^solving a system of nonlinear equations (parallel Bratu example);
16: Concepts: DA^using distributed arrays;
17: Processors: n
18: T*/
20: /* ------------------------------------------------------------------------
22: Modified Solid Fuel Ignition problem. This problem is modeled by
23: the partial differential equation
25: -Laplacian u - kappa*\PartDer{u}{x} - lambda*exp(u) = 0,
27: where
29: 0 < x,y < 1,
30:
31: with boundary conditions
32:
33: u = 0 for x = 0, x = 1, y = 0, y = 1.
34:
35: A finite difference approximation with the usual 5-point stencil
36: is used to discretize the boundary value problem to obtain a nonlinear
37: system of equations.
39: ------------------------------------------------------------------------- */
41: /*
42: Include "petscda.h" so that we can use distributed arrays (DAs).
43: Include "petscsnes.h" so that we can use SNES solvers. Note that this
44: file automatically includes:
45: petsc.h - base PETSc routines petscvec.h - vectors
46: petscsys.h - system routines petscmat.h - matrices
47: petscis.h - index sets petscksp.h - Krylov subspace methods
48: petscviewer.h - viewers petscpc.h - preconditioners
49: petscsles.h - linear solvers
50: */
51: #include petscda.h
52: #include petscsnes.h
54: /*
55: User-defined application context - contains data needed by the
56: application-provided call-back routines, FormJacobian() and
57: FormFunction().
58: */
59: typedef struct {
60: PetscReal param; /* test problem parameter */
61: PetscReal param2; /* test problem parameter */
62: int mx,my; /* discretization in x, y directions */
63: Vec localX,localF; /* ghosted local vector */
64: DA da; /* distributed array data structure */
65: int rank; /* processor rank */
66: } AppCtx;
68: /*
69: User-defined routines
70: */
71: extern int FormFunction(SNES,Vec,Vec,void*),FormInitialGuess(AppCtx*,Vec);
72: extern int FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
76: int main(int argc,char **argv)
77: {
78: SNES snes; /* nonlinear solver */
79: Vec x,r; /* solution, residual vectors */
80: Mat J; /* Jacobian matrix */
81: AppCtx user; /* user-defined work context */
82: int its; /* iterations for convergence */
83: int Nx,Ny; /* number of preocessors in x- and y- directions */
84: PetscTruth matrix_free; /* flag - 1 indicates matrix-free version */
85: int size; /* number of processors */
86: int m,N,ierr;
87: PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
88: PetscReal bratu_kappa_max = 10000,bratu_kappa_min = 0.;
90: PetscInitialize(&argc,&argv,(char *)0,help);
91: MPI_Comm_rank(PETSC_COMM_WORLD,&user.rank);
93: /*
94: Initialize problem parameters
95: */
96: user.mx = 4; user.my = 4; user.param = 6.0; user.param2 = 0.0;
97: PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
98: PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
99: PetscOptionsGetReal(PETSC_NULL,"-lambda",&user.param,PETSC_NULL);
100: if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
101: SETERRQ(1,"Lambda is out of range");
102: }
103: PetscOptionsGetReal(PETSC_NULL,"-kappa",&user.param2,PETSC_NULL);
104: if (user.param2 >= bratu_kappa_max || user.param2 < bratu_kappa_min) {
105: SETERRQ(1,"Kappa is out of range");
106: }
107: PetscPrintf(PETSC_COMM_WORLD,"Solving the Bratu problem with lambda=%g, kappa=%g\n",user.param,user.param2);
109: N = user.mx*user.my;
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Create nonlinear solver context
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: SNESCreate(PETSC_COMM_WORLD,&snes);
117: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118: Create vector data structures; set function evaluation routine
119: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
121: /*
122: Create distributed array (DA) to manage parallel grid and vectors
123: */
124: MPI_Comm_size(PETSC_COMM_WORLD,&size);
125: Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
126: PetscOptionsGetInt(PETSC_NULL,"-Nx",&Nx,PETSC_NULL);
127: PetscOptionsGetInt(PETSC_NULL,"-Ny",&Ny,PETSC_NULL);
128: if (Nx*Ny != size && (Nx != PETSC_DECIDE || Ny != PETSC_DECIDE))
129: SETERRQ(1,"Incompatible number of processors: Nx * Ny != size");
130: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,user.mx,user.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&user.da);
132: /*
133: Visualize the distribution of the array across the processors
134: */
135: /* DAView(user.da,PETSC_VIEWER_DRAW_WORLD); */
138: /*
139: Extract global and local vectors from DA; then duplicate for remaining
140: vectors that are the same types
141: */
142: DACreateGlobalVector(user.da,&x);
143: DACreateLocalVector(user.da,&user.localX);
144: VecDuplicate(x,&r);
145: VecDuplicate(user.localX,&user.localF);
147: /*
148: Set function evaluation routine and vector
149: */
150: SNESSetFunction(snes,r,FormFunction,(void*)&user);
152: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153: Create matrix data structure; set Jacobian evaluation routine
154: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156: /*
157: Set Jacobian matrix data structure and default Jacobian evaluation
158: routine. User can override with:
159: -snes_fd : default finite differencing approximation of Jacobian
160: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
161: (unless user explicitly sets preconditioner)
162: -snes_mf_operator : form preconditioning matrix as set by the user,
163: but use matrix-free approx for Jacobian-vector
164: products within Newton-Krylov method
166: Note: For the parallel case, vectors and matrices MUST be partitioned
167: accordingly. When using distributed arrays (DAs) to create vectors,
168: the DAs determine the problem partitioning. We must explicitly
169: specify the local matrix dimensions upon its creation for compatibility
170: with the vector distribution. Thus, the generic MatCreate() routine
171: is NOT sufficient when working with distributed arrays.
173: Note: Here we only approximately preallocate storage space for the
174: Jacobian. See the users manual for a discussion of better techniques
175: for preallocating matrix memory.
176: */
177: PetscOptionsHasName(PETSC_NULL,"-snes_mf",&matrix_free);
178: if (!matrix_free) {
179: if (size == 1) {
180: MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,PETSC_NULL,&J);
181: } else {
182: VecGetLocalSize(x,&m);
183: MatCreateMPIAIJ(PETSC_COMM_WORLD,m,m,N,N,5,PETSC_NULL,3,PETSC_NULL,&J);
184: }
185: SNESSetJacobian(snes,J,J,FormJacobian,&user);
186: }
188: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189: Customize nonlinear solver; set runtime options
190: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192: /*
193: Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
194: */
195: SNESSetFromOptions(snes);
197: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198: Evaluate initial guess; then solve nonlinear system
199: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200: /*
201: Note: The user should initialize the vector, x, with the initial guess
202: for the nonlinear solver prior to calling SNESSolve(). In particular,
203: to employ an initial guess of zero, the user should explicitly set
204: this vector to zero by calling VecSet().
205: */
206: FormInitialGuess(&user,x);
207: SNESSolve(snes,x,&its);
208: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %d\n",its);
210: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: Free work space. All PETSc objects should be destroyed when they
212: are no longer needed.
213: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
215: if (!matrix_free) {
216: MatDestroy(J);
217: }
218: VecDestroy(user.localX); VecDestroy(x);
219: VecDestroy(user.localF); VecDestroy(r);
220: SNESDestroy(snes); DADestroy(user.da);
221: PetscFinalize();
223: return 0;
224: }
225: /* ------------------------------------------------------------------- */
228: /*
229: FormInitialGuess - Forms initial approximation.
231: Input Parameters:
232: user - user-defined application context
233: X - vector
235: Output Parameter:
236: X - vector
237: */
238: int FormInitialGuess(AppCtx *user,Vec X)
239: {
240: int i,j,row,mx,my,ierr,xs,ys,xm,ym,gxm,gym,gxs,gys;
241: PetscReal one = 1.0,lambda,temp1,temp,hx,hy,hxdhy,hydhx,sc;
242: PetscScalar *x;
243: Vec localX = user->localX;
245: mx = user->mx; my = user->my; lambda = user->param;
246: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
247: sc = hx*hy*lambda; hxdhy = hx/hy; hydhx = hy/hx;
248: temp1 = lambda/(lambda + one);
250: /*
251: Get a pointer to vector data.
252: - For default PETSc vectors,VecGetArray() returns a pointer to
253: the data array. Otherwise, the routine is implementation dependent.
254: - You MUST call VecRestoreArray() when you no longer need access to
255: the array.
256: */
257: VecGetArray(localX,&x);
259: /*
260: Get local grid boundaries (for 2-dimensional DA):
261: xs, ys - starting grid indices (no ghost points)
262: xm, ym - widths of local grid (no ghost points)
263: gxs, gys - starting grid indices (including ghost points)
264: gxm, gym - widths of local grid (including ghost points)
265: */
266: DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
267: DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL);
269: /*
270: Compute initial guess over the locally owned part of the grid
271: */
272: for (j=ys; j<ys+ym; j++) {
273: temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
274: for (i=xs; i<xs+xm; i++) {
275: row = i - gxs + (j - gys)*gxm;
276: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
277: x[row] = 0.0;
278: continue;
279: }
280: x[row] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
281: }
282: }
284: /*
285: Restore vector
286: */
287: VecRestoreArray(localX,&x);
289: /*
290: Insert values into global vector
291: */
292: DALocalToGlobal(user->da,localX,INSERT_VALUES,X);
293: return 0;
294: }
295: /* ------------------------------------------------------------------- */
298: /*
299: FormFunction - Evaluates nonlinear function, F(x).
301: Input Parameters:
302: . snes - the SNES context
303: . X - input vector
304: . ptr - optional user-defined context, as set by SNESSetFunction()
306: Output Parameter:
307: . F - function vector
308: */
309: int FormFunction(SNES snes,Vec X,Vec F,void *ptr)
310: {
311: AppCtx *user = (AppCtx*)ptr;
312: int ierr,i,j,row,mx,my,xs,ys,xm,ym,gxs,gys,gxm,gym;
313: PetscReal two = 2.0,one = 1.0,half = 0.5;
314: PetscReal lambda,hx,hy,hxdhy,hydhx,sc;
315: PetscScalar u,ux,uxx,uyy,*x,*f,kappa;
316: Vec localX = user->localX,localF = user->localF;
318: mx = user->mx; my = user->my; lambda = user->param;
319: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
320: sc = hx*hy*lambda; hxdhy = hx/hy; hydhx = hy/hx;
321: kappa = user->param2;
323: /*
324: Scatter ghost points to local vector, using the 2-step process
325: DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
326: By placing code between these two statements, computations can be
327: done while messages are in transition.
328: */
329: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
330: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
332: /*
333: Get pointers to vector data
334: */
335: VecGetArray(localX,&x);
336: VecGetArray(localF,&f);
338: /*
339: Get local grid boundaries
340: */
341: DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
342: DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL);
344: /*
345: Compute function over the locally owned part of the grid
346: */
347: for (j=ys; j<ys+ym; j++) {
348: row = (j - gys)*gxm + xs - gxs - 1;
349: for (i=xs; i<xs+xm; i++) {
350: row++;
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: ux = (x[row+1] - x[row-1])*half*hy;
357: uxx = (two*u - x[row-1] - x[row+1])*hydhx;
358: uyy = (two*u - x[row-gxm] - x[row+gxm])*hxdhy;
359: f[row] = uxx + uyy - kappa*ux - sc*exp(u);
360: }
361: }
363: /*
364: Restore vectors
365: */
366: VecRestoreArray(localX,&x);
367: VecRestoreArray(localF,&f);
369: /*
370: Insert values into global vector
371: */
372: DALocalToGlobal(user->da,localF,INSERT_VALUES,F);
373: PetscLogFlops(11*ym*xm);
374: return 0;
375: }
376: /* ------------------------------------------------------------------- */
379: /*
380: FormJacobian - Evaluates Jacobian matrix.
382: Input Parameters:
383: . snes - the SNES context
384: . x - input vector
385: . ptr - optional user-defined context, as set by SNESSetJacobian()
387: Output Parameters:
388: . A - Jacobian matrix
389: . B - optionally different preconditioning matrix
390: . flag - flag indicating matrix structure
392: Notes:
393: Due to grid point reordering with DAs, we must always work
394: with the local grid points, and then transform them to the new
395: global numbering with the "ltog" mapping (via DAGetGlobalIndices()).
396: We cannot work directly with the global numbers for the original
397: uniprocessor grid!
398: */
399: int FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
400: {
401: AppCtx *user = (AppCtx*)ptr; /* user-defined application context */
402: Mat jac = *B; /* Jacobian matrix */
403: Vec localX = user->localX; /* local vector */
404: int *ltog; /* local-to-global mapping */
405: int ierr,i,j,row,mx,my,col[5];
406: int nloc,xs,ys,xm,ym,gxs,gys,gxm,gym,grow;
407: PetscScalar two = 2.0,one = 1.0,lambda,v[5],hx,hy,hxdhy,hydhx,sc,*x;
409: mx = user->mx; my = user->my; lambda = user->param;
410: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
411: sc = hx*hy; hxdhy = hx/hy; hydhx = hy/hx;
413: /*
414: Scatter ghost points to local vector,using the 2-step process
415: DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
416: By placing code between these two statements, computations can be
417: done while messages are in transition.
418: */
419: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
420: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
422: /*
423: Get pointer to vector data
424: */
425: VecGetArray(localX,&x);
427: /*
428: Get local grid boundaries
429: */
430: DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
431: DAGetGhostCorners(user->da,&gxs,&gys,PETSC_NULL,&gxm,&gym,PETSC_NULL);
433: /*
434: Get the global node numbers for all local nodes, including ghost points
435: */
436: DAGetGlobalIndices(user->da,&nloc,<og);
438: /*
439: Compute entries for the locally owned part of the Jacobian.
440: - Currently, all PETSc parallel matrix formats are partitioned by
441: contiguous chunks of rows across the processors. The "grow"
442: parameter computed below specifies the global row number
443: corresponding to each local grid point.
444: - Each processor needs to insert only elements that it owns
445: locally (but any non-local elements will be sent to the
446: appropriate processor during matrix assembly).
447: - Always specify global row and columns of matrix entries.
448: - Here, we set all entries for a particular row at once.
449: */
450: for (j=ys; j<ys+ym; j++) {
451: row = (j - gys)*gxm + xs - gxs - 1;
452: for (i=xs; i<xs+xm; i++) {
453: row++;
454: grow = ltog[row];
455: /* boundary points */
456: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
457: MatSetValues(jac,1,&grow,1,&grow,&one,INSERT_VALUES);
458: continue;
459: }
460: /* interior grid points */
461: v[0] = -hxdhy; col[0] = ltog[row - gxm];
462: v[1] = -hydhx; col[1] = ltog[row - 1];
463: v[2] = two*(hydhx + hxdhy) - sc*lambda*exp(x[row]); col[2] = grow;
464: v[3] = -hydhx; col[3] = ltog[row + 1];
465: v[4] = -hxdhy; col[4] = ltog[row + gxm];
466: MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
467: }
468: }
470: /*
471: Assemble matrix, using the 2-step process:
472: MatAssemblyBegin(), MatAssemblyEnd().
473: By placing code between these two statements, computations can be
474: done while messages are in transition.
475: */
476: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
477: VecRestoreArray(localX,&x);
478: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
480: /*
481: Set flag to indicate that the Jacobian matrix retains an identical
482: nonzero structure throughout all nonlinear iterations (although the
483: values of the entries change). Thus, we can save some work in setting
484: up the preconditioner (e.g., no need to redo symbolic factorization for
485: ILU/ICC preconditioners).
486: - If the nonzero structure of the matrix is different during
487: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
488: must be used instead. If you are unsure whether the matrix
489: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
490: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
491: believes your assertion and does not check the structure
492: of the matrix. If you erroneously claim that the structure
493: is the same when it actually is not, the new preconditioner
494: will not function correctly. Thus, use this optimization
495: feature with caution!
496: */
497: *flag = SAME_NONZERO_PATTERN;
498: /*
499: Tell the matrix we will never add a new nonzero location to the
500: matrix. If we do it will generate an error.
501: */
502: /* MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR); */
503: return 0;
504: }