Actual source code: ex14.c
1: /*$Id: ex14.c,v 1.23 2001/08/07 21:31:17 bsmith Exp $*/
3: /* Program usage: mpirun -np <procs> ex14 [-help] [all PETSc options] */
5: static char help[] = "Bratu nonlinear PDE in 3d.n
6: We solve the Bratu (SFI - solid fuel ignition) problem in a 3D rectangularn
7: domain, using distributed arrays (DAs) to partition the parallel grid.n
8: The command line options include:n
9: -par <parameter>, where <parameter> indicates the problem's nonlinearityn
10: problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)nn";
12: /*T
13: Concepts: SNES^parallel Bratu example
14: Concepts: DA^using distributed arrays;
15: Processors: n
16: T*/
18: /* ------------------------------------------------------------------------
20: Solid Fuel Ignition (SFI) problem. This problem is modeled by
21: the partial differential equation
22:
23: -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
24:
25: with boundary conditions
26:
27: u = 0 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1
28:
29: A finite difference approximation with the usual 7-point stencil
30: is used to discretize the boundary value problem to obtain a nonlinear
31: system of equations.
34: ------------------------------------------------------------------------- */
36: /*
37: Include "petscda.h" so that we can use distributed arrays (DAs).
38: Include "petscsnes.h" so that we can use SNES solvers. Note that this
39: file automatically includes:
40: petsc.h - base PETSc routines petscvec.h - vectors
41: petscsys.h - system routines petscmat.h - matrices
42: petscis.h - index sets petscksp.h - Krylov subspace methods
43: petscviewer.h - viewers petscpc.h - preconditioners
44: petscsles.h - linear solvers
45: */
46: #include petscda.h
47: #include petscsnes.h
50: /*
51: User-defined application context - contains data needed by the
52: application-provided call-back routines, FormJacobian() and
53: FormFunction().
54: */
55: typedef struct {
56: PetscReal param; /* test problem parameter */
57: DA da; /* distributed array data structure */
58: } AppCtx;
60: /*
61: User-defined routines
62: */
63: extern int FormFunction(SNES,Vec,Vec,void*),FormInitialGuess(AppCtx*,Vec);
64: extern int FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
66: int main(int argc,char **argv)
67: {
68: SNES snes; /* nonlinear solver */
69: Vec x,r; /* solution, residual vectors */
70: Mat J; /* Jacobian matrix */
71: AppCtx user; /* user-defined work context */
72: int its; /* iterations for convergence */
73: PetscTruth matrix_free,coloring;
74: int ierr;
75: PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.,fnorm;
76: MatFDColoring matfdcoloring;
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Initialize program
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: PetscInitialize(&argc,&argv,(char *)0,help);
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Initialize problem parameters
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: user.param = 6.0;
88: PetscOptionsGetReal(PETSC_NULL,"-par",&user.param,PETSC_NULL);
89: if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
90: SETERRQ(1,"Lambda is out of range");
91: }
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Create nonlinear solver context
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: SNESCreate(PETSC_COMM_WORLD,&snes);
98: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99: Create distributed array (DA) to manage parallel grid and vectors
100: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101: DACreate3d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,4,4,4,PETSC_DECIDE,PETSC_DECIDE,
102: PETSC_DECIDE,1,1,PETSC_NULL,PETSC_NULL,PETSC_NULL,&user.da);
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Extract global vectors from DA; then duplicate for remaining
106: vectors that are the same types
107: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: DACreateGlobalVector(user.da,&x);
109: VecDuplicate(x,&r);
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Set function evaluation routine and vector
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: SNESSetFunction(snes,r,FormFunction,(void*)&user);
116: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: Create matrix data structure; set Jacobian evaluation routine
119: Set Jacobian matrix data structure and default Jacobian evaluation
120: routine. User can override with:
121: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
122: (unless user explicitly sets preconditioner)
123: -snes_mf_operator : form preconditioning matrix as set by the user,
124: but use matrix-free approx for Jacobian-vector
125: products within Newton-Krylov method
126: -fdcoloring : using finite differences with coloring to compute the Jacobian
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129: PetscOptionsHasName(PETSC_NULL,"-snes_mf",&matrix_free);
130: PetscOptionsHasName(PETSC_NULL,"-fdcoloring",&coloring);
131: if (!matrix_free) {
132: if (coloring) {
133: ISColoring iscoloring;
135: DAGetColoring(user.da,IS_COLORING_LOCAL,&iscoloring);
136: DAGetMatrix(user.da,MATMPIAIJ,&J);
137: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
138: ISColoringDestroy(iscoloring);
139: MatFDColoringSetFunction(matfdcoloring,(int (*)(void))FormFunction,&user);
140: MatFDColoringSetFromOptions(matfdcoloring);
141: SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,matfdcoloring);
142: } else {
143: DAGetMatrix(user.da,MATMPIAIJ,&J);
144: SNESSetJacobian(snes,J,J,FormJacobian,&user);
145: }
146: }
148: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149: Customize nonlinear solver; set runtime options
150: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151: SNESSetFromOptions(snes);
153: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154: Evaluate initial guess
155: Note: The user should initialize the vector, x, with the initial guess
156: for the nonlinear solver prior to calling SNESSolve(). In particular,
157: to employ an initial guess of zero, the user should explicitly set
158: this vector to zero by calling VecSet().
159: */
160: FormInitialGuess(&user,x);
162: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163: Solve nonlinear system
164: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165: SNESSolve(snes,x,&its);
167: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168: Explicitly check norm of the residual of the solution
169: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170: FormFunction(snes,x,r,(void *)&user);
171: VecNorm(r,NORM_2,&fnorm);
172: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %d fnorm %gn",its,fnorm);
174: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175: Free work space. All PETSc objects should be destroyed when they
176: are no longer needed.
177: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
179: if (!matrix_free) {
180: MatDestroy(J);
181: }
182: if (coloring) {
183: MatFDColoringDestroy(matfdcoloring);
184: }
185: VecDestroy(x);
186: VecDestroy(r);
187: SNESDestroy(snes);
188: DADestroy(user.da);
189: PetscFinalize();
191: return(0);
192: }
193: /* ------------------------------------------------------------------- */
194: /*
195: FormInitialGuess - Forms initial approximation.
197: Input Parameters:
198: user - user-defined application context
199: X - vector
201: Output Parameter:
202: X - vector
203: */
204: int FormInitialGuess(AppCtx *user,Vec X)
205: {
206: int i,j,k,Mx,My,Mz,ierr,xs,ys,zs,xm,ym,zm;
207: PetscReal lambda,temp1,hx,hy,hz,tempk,tempj;
208: PetscScalar ***x;
211: DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
212: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
214: lambda = user->param;
215: hx = 1.0/(PetscReal)(Mx-1);
216: hy = 1.0/(PetscReal)(My-1);
217: hz = 1.0/(PetscReal)(Mz-1);
218: temp1 = lambda/(lambda + 1.0);
220: /*
221: Get a pointer to vector data.
222: - For default PETSc vectors, VecGetArray() returns a pointer to
223: the data array. Otherwise, the routine is implementation dependent.
224: - You MUST call VecRestoreArray() when you no longer need access to
225: the array.
226: */
227: DAVecGetArray(user->da,X,(void**)&x);
229: /*
230: Get local grid boundaries (for 3-dimensional DA):
231: xs, ys, zs - starting grid indices (no ghost points)
232: xm, ym, zm - widths of local grid (no ghost points)
234: */
235: DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);
237: /*
238: Compute initial guess over the locally owned part of the grid
239: */
240: for (k=zs; k<zs+zm; k++) {
241: tempk = (PetscReal)(PetscMin(k,Mz-k-1))*hz;
242: for (j=ys; j<ys+ym; j++) {
243: tempj = PetscMin((PetscReal)(PetscMin(j,My-j-1))*hy,tempk);
244: for (i=xs; i<xs+xm; i++) {
245: if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
246: /* boundary conditions are all zero Dirichlet */
247: x[k][j][i] = 0.0;
248: } else {
249: x[k][j][i] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,tempj));
250: }
251: }
252: }
253: }
255: /*
256: Restore vector
257: */
258: DAVecRestoreArray(user->da,X,(void**)&x);
259: return(0);
260: }
261: /* ------------------------------------------------------------------- */
262: /*
263: FormFunction - Evaluates nonlinear function, F(x).
265: Input Parameters:
266: . snes - the SNES context
267: . X - input vector
268: . ptr - optional user-defined context, as set by SNESSetFunction()
270: Output Parameter:
271: . F - function vector
272: */
273: int FormFunction(SNES snes,Vec X,Vec F,void *ptr)
274: {
275: AppCtx *user = (AppCtx*)ptr;
276: int ierr,i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
277: PetscReal two = 2.0,lambda,hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc;
278: PetscScalar u_north,u_south,u_east,u_west,u_up,u_down,u,u_xx,u_yy,u_zz,***x,***f;
279: Vec localX;
282: DAGetLocalVector(user->da,&localX);
283: DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
284: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
286: lambda = user->param;
287: hx = 1.0/(PetscReal)(Mx-1);
288: hy = 1.0/(PetscReal)(My-1);
289: hz = 1.0/(PetscReal)(Mz-1);
290: sc = hx*hy*hz*lambda;
291: hxhzdhy = hx*hz/hy;
292: hyhzdhx = hy*hz/hx;
293: hxhydhz = hx*hy/hz;
295: /*
296: Scatter ghost points to local vector,using the 2-step process
297: DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
298: By placing code between these two statements, computations can be
299: done while messages are in transition.
300: */
301: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
302: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
304: /*
305: Get pointers to vector data
306: */
307: DAVecGetArray(user->da,localX,(void**)&x);
308: DAVecGetArray(user->da,F,(void**)&f);
310: /*
311: Get local grid boundaries
312: */
313: DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);
315: /*
316: Compute function over the locally owned part of the grid
317: */
318: for (k=zs; k<zs+zm; k++) {
319: for (j=ys; j<ys+ym; j++) {
320: for (i=xs; i<xs+xm; i++) {
321: if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
322: f[k][j][i] = x[k][j][i];
323: } else {
324: u = x[k][j][i];
325: u_east = x[k][j][i+1];
326: u_west = x[k][j][i-1];
327: u_north = x[k][j+1][i];
328: u_south = x[k][j-1][i];
329: u_up = x[k+1][j][i];
330: u_down = x[k-1][j][i];
331: u_xx = (-u_east + two*u - u_west)*hyhzdhx;
332: u_yy = (-u_north + two*u - u_south)*hxhzdhy;
333: u_zz = (-u_up + two*u - u_down)*hxhydhz;
334: f[k][j][i] = u_xx + u_yy + u_zz - sc*PetscExpScalar(u);
335: }
336: }
337: }
338: }
340: /*
341: Restore vectors
342: */
343: DAVecRestoreArray(user->da,localX,(void**)&x);
344: DAVecRestoreArray(user->da,F,(void**)&f);
345: DARestoreLocalVector(user->da,&localX);
346: PetscLogFlops(11*ym*xm);
347: return(0);
348: }
349: /* ------------------------------------------------------------------- */
350: /*
351: FormJacobian - Evaluates Jacobian matrix.
353: Input Parameters:
354: . snes - the SNES context
355: . x - input vector
356: . ptr - optional user-defined context, as set by SNESSetJacobian()
358: Output Parameters:
359: . A - Jacobian matrix
360: . B - optionally different preconditioning matrix
361: . flag - flag indicating matrix structure
363: */
364: int FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
365: {
366: AppCtx *user = (AppCtx*)ptr; /* user-defined application context */
367: Mat jac = *B; /* Jacobian matrix */
368: Vec localX;
369: int ierr,i,j,k,Mx,My,Mz;
370: MatStencil col[7],row;
371: int xs,ys,zs,xm,ym,zm;
372: PetscScalar lambda,v[7],hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc,***x;
376: DAGetLocalVector(user->da,&localX);
377: DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
378: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
380: lambda = user->param;
381: hx = 1.0/(PetscReal)(Mx-1);
382: hy = 1.0/(PetscReal)(My-1);
383: hz = 1.0/(PetscReal)(Mz-1);
384: sc = hx*hy*hz*lambda;
385: hxhzdhy = hx*hz/hy;
386: hyhzdhx = hy*hz/hx;
387: hxhydhz = hx*hy/hz;
389: /*
390: Scatter ghost points to local vector, using the 2-step process
391: DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
392: By placing code between these two statements, computations can be
393: done while messages are in transition.
394: */
395: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
396: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
398: /*
399: Get pointer to vector data
400: */
401: DAVecGetArray(user->da,localX,(void**)&x);
403: /*
404: Get local grid boundaries
405: */
406: DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);
408: /*
409: Compute entries for the locally owned part of the Jacobian.
410: - Currently, all PETSc parallel matrix formats are partitioned by
411: contiguous chunks of rows across the processors.
412: - Each processor needs to insert only elements that it owns
413: locally (but any non-local elements will be sent to the
414: appropriate processor during matrix assembly).
415: - Here, we set all entries for a particular row at once.
416: - We can set matrix entries either using either
417: MatSetValuesLocal() or MatSetValues(), as discussed above.
418: */
419: for (k=zs; k<zs+zm; k++) {
420: for (j=ys; j<ys+ym; j++) {
421: for (i=xs; i<xs+xm; i++) {
422: row.k = k; row.j = j; row.i = i;
423: /* boundary points */
424: if (i == 0 || j == 0 || k == 0|| i == Mx-1 || j == My-1 || k == Mz-1) {
425: v[0] = 1.0;
426: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
427: } else {
428: /* interior grid points */
429: v[0] = -hxhydhz; col[0].k=k-1;col[0].j=j; col[0].i = i;
430: v[1] = -hxhzdhy; col[1].k=k; col[1].j=j-1;col[1].i = i;
431: v[2] = -hyhzdhx; col[2].k=k; col[2].j=j; col[2].i = i-1;
432: v[3] = 2.0*(hyhzdhx+hxhzdhy+hxhydhz)-sc*PetscExpScalar(x[k][j][i]);col[3].k=row.k;col[3].j=row.j;col[3].i = row.i;
433: v[4] = -hyhzdhx; col[4].k=k; col[4].j=j; col[4].i = i+1;
434: v[5] = -hxhzdhy; col[5].k=k; col[5].j=j+1;col[5].i = i;
435: v[6] = -hxhydhz; col[6].k=k+1;col[6].j=j; col[6].i = i;
436: MatSetValuesStencil(jac,1,&row,7,col,v,INSERT_VALUES);
437: }
438: }
439: }
440: }
441: DAVecRestoreArray(user->da,localX,(void**)&x);
442: DARestoreLocalVector(user->da,&localX);
444: /*
445: Assemble matrix, using the 2-step process:
446: MatAssemblyBegin(), MatAssemblyEnd().
447: */
448: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
449: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
451: /*
452: Normally since the matrix has already been assembled above; this
453: would do nothing. But in the matrix free mode -snes_mf_operator
454: this tells the "matrix-free" matrix that a new linear system solve
455: is about to be done.
456: */
458: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
459: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
461: /*
462: Set flag to indicate that the Jacobian matrix retains an identical
463: nonzero structure throughout all nonlinear iterations (although the
464: values of the entries change). Thus, we can save some work in setting
465: up the preconditioner (e.g., no need to redo symbolic factorization for
466: ILU/ICC preconditioners).
467: - If the nonzero structure of the matrix is different during
468: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
469: must be used instead. If you are unsure whether the matrix
470: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
471: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
472: believes your assertion and does not check the structure
473: of the matrix. If you erroneously claim that the structure
474: is the same when it actually is not, the new preconditioner
475: will not function correctly. Thus, use this optimization
476: feature with caution!
477: */
478: *flag = SAME_NONZERO_PATTERN;
481: /*
482: Tell the matrix we will never add a new nonzero location to the
483: matrix. If we do, it will generate an error.
484: */
485: MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR);
486: return(0);
487: }