Actual source code: ex14.c
petsc-3.4.0 2013-05-13
2: static char help[] = "Bratu nonlinear PDE in 3d.\n\
3: We solve the Bratu (SFI - solid fuel ignition) problem in a 3D rectangular\n\
4: domain, using distributed arrays (DMDAs) to partition the parallel grid.\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\n";
9: /*T
10: Concepts: SNES^parallel Bratu example
11: Concepts: DMDA^using distributed arrays;
12: Processors: n
13: T*/
15: /* ------------------------------------------------------------------------
17: Solid Fuel Ignition (SFI) problem. This problem is modeled by
18: the partial differential equation
20: -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
22: with boundary conditions
24: u = 0 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1
26: A finite difference approximation with the usual 7-point stencil
27: is used to discretize the boundary value problem to obtain a nonlinear
28: system of equations.
31: ------------------------------------------------------------------------- */
33: /*
34: Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
35: Include "petscsnes.h" so that we can use SNES solvers. Note that this
36: file automatically includes:
37: petscsys.h - base PETSc routines petscvec.h - vectors
38: petscmat.h - matrices
39: petscis.h - index sets petscksp.h - Krylov subspace methods
40: petscviewer.h - viewers petscpc.h - preconditioners
41: petscksp.h - linear solvers
42: */
43: #include <petscdmda.h>
44: #include <petscsnes.h>
47: /*
48: User-defined application context - contains data needed by the
49: application-provided call-back routines, FormJacobian() and
50: FormFunction().
51: */
52: typedef struct {
53: PetscReal param; /* test problem parameter */
54: DM da; /* distributed array data structure */
55: } AppCtx;
57: /*
58: User-defined routines
59: */
60: extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*),FormInitialGuess(AppCtx*,Vec);
61: extern PetscErrorCode FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
65: int main(int argc,char **argv)
66: {
67: SNES snes; /* nonlinear solver */
68: Vec x,r; /* solution, residual vectors */
69: Mat J; /* Jacobian matrix */
70: AppCtx user; /* user-defined work context */
71: PetscInt its; /* iterations for convergence */
72: PetscBool matrix_free = PETSC_FALSE,coloring = PETSC_FALSE;
74: PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.,fnorm;
76: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77: Initialize program
78: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80: PetscInitialize(&argc,&argv,(char*)0,help);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Initialize problem parameters
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: user.param = 6.0;
86: PetscOptionsGetReal(NULL,"-par",&user.param,NULL);
87: if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) SETERRQ(PETSC_COMM_SELF,1,"Lambda is out of range");
89: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90: Create nonlinear solver context
91: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: SNESCreate(PETSC_COMM_WORLD,&snes);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Create distributed array (DMDA) to manage parallel grid and vectors
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: DMDACreate3d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,-4,PETSC_DECIDE,PETSC_DECIDE,
98: PETSC_DECIDE,1,1,NULL,NULL,NULL,&user.da);
100: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: Extract global vectors from DMDA; then duplicate for remaining
102: vectors that are the same types
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: DMCreateGlobalVector(user.da,&x);
105: VecDuplicate(x,&r);
107: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108: Set function evaluation routine and vector
109: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110: SNESSetFunction(snes,r,FormFunction,(void*)&user);
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Create matrix data structure; set Jacobian evaluation routine
115: Set Jacobian matrix data structure and default Jacobian evaluation
116: routine. User can override with:
117: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
118: (unless user explicitly sets preconditioner)
119: -snes_mf_operator : form preconditioning matrix as set by the user,
120: but use matrix-free approx for Jacobian-vector
121: products within Newton-Krylov method
122: -fdcoloring : using finite differences with coloring to compute the Jacobian
124: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125: PetscOptionsGetBool(NULL,"-snes_mf",&matrix_free,NULL);
126: PetscOptionsGetBool(NULL,"-fdcoloring",&coloring,NULL);
127: if (!matrix_free) {
128: DMCreateMatrix(user.da,MATAIJ,&J);
129: if (coloring) {
130: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,0);
131: } else {
132: SNESSetJacobian(snes,J,J,FormJacobian,&user);
133: }
134: }
136: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137: Customize nonlinear solver; set runtime options
138: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139: SNESSetFromOptions(snes);
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Evaluate initial guess
143: Note: The user should initialize the vector, x, with the initial guess
144: for the nonlinear solver prior to calling SNESSolve(). In particular,
145: to employ an initial guess of zero, the user should explicitly set
146: this vector to zero by calling VecSet().
147: */
148: FormInitialGuess(&user,x);
150: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151: Solve nonlinear system
152: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153: SNESSolve(snes,NULL,x);
154: SNESGetIterationNumber(snes,&its);
156: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157: Explicitly check norm of the residual of the solution
158: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159: FormFunction(snes,x,r,(void*)&user);
160: VecNorm(r,NORM_2,&fnorm);
161: PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D fnorm %G\n",its,fnorm);
163: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164: Free work space. All PETSc objects should be destroyed when they
165: are no longer needed.
166: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168: if (!matrix_free) {
169: MatDestroy(&J);
170: }
171: VecDestroy(&x);
172: VecDestroy(&r);
173: SNESDestroy(&snes);
174: DMDestroy(&user.da);
175: PetscFinalize();
176: return(0);
177: }
178: /* ------------------------------------------------------------------- */
181: /*
182: FormInitialGuess - Forms initial approximation.
184: Input Parameters:
185: user - user-defined application context
186: X - vector
188: Output Parameter:
189: X - vector
190: */
191: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
192: {
193: PetscInt i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
195: PetscReal lambda,temp1,hx,hy,hz,tempk,tempj;
196: PetscScalar ***x;
199: DMDAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
201: lambda = user->param;
202: hx = 1.0/(PetscReal)(Mx-1);
203: hy = 1.0/(PetscReal)(My-1);
204: hz = 1.0/(PetscReal)(Mz-1);
205: temp1 = lambda/(lambda + 1.0);
207: /*
208: Get a pointer to vector data.
209: - For default PETSc vectors, VecGetArray() returns a pointer to
210: the data array. Otherwise, the routine is implementation dependent.
211: - You MUST call VecRestoreArray() when you no longer need access to
212: the array.
213: */
214: DMDAVecGetArray(user->da,X,&x);
216: /*
217: Get local grid boundaries (for 3-dimensional DMDA):
218: xs, ys, zs - starting grid indices (no ghost points)
219: xm, ym, zm - widths of local grid (no ghost points)
221: */
222: DMDAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);
224: /*
225: Compute initial guess over the locally owned part of the grid
226: */
227: for (k=zs; k<zs+zm; k++) {
228: tempk = (PetscReal)(PetscMin(k,Mz-k-1))*hz;
229: for (j=ys; j<ys+ym; j++) {
230: tempj = PetscMin((PetscReal)(PetscMin(j,My-j-1))*hy,tempk);
231: for (i=xs; i<xs+xm; i++) {
232: if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
233: /* boundary conditions are all zero Dirichlet */
234: x[k][j][i] = 0.0;
235: } else {
236: x[k][j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,tempj));
237: }
238: }
239: }
240: }
242: /*
243: Restore vector
244: */
245: DMDAVecRestoreArray(user->da,X,&x);
246: return(0);
247: }
248: /* ------------------------------------------------------------------- */
251: /*
252: FormFunction - Evaluates nonlinear function, F(x).
254: Input Parameters:
255: . snes - the SNES context
256: . X - input vector
257: . ptr - optional user-defined context, as set by SNESSetFunction()
259: Output Parameter:
260: . F - function vector
261: */
262: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
263: {
264: AppCtx *user = (AppCtx*)ptr;
266: PetscInt i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
267: PetscReal two = 2.0,lambda,hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc;
268: PetscScalar u_north,u_south,u_east,u_west,u_up,u_down,u,u_xx,u_yy,u_zz,***x,***f;
269: Vec localX;
272: DMGetLocalVector(user->da,&localX);
273: DMDAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
274: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
276: lambda = user->param;
277: hx = 1.0/(PetscReal)(Mx-1);
278: hy = 1.0/(PetscReal)(My-1);
279: hz = 1.0/(PetscReal)(Mz-1);
280: sc = hx*hy*hz*lambda;
281: hxhzdhy = hx*hz/hy;
282: hyhzdhx = hy*hz/hx;
283: hxhydhz = hx*hy/hz;
285: /*
286: Scatter ghost points to local vector,using the 2-step process
287: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
288: By placing code between these two statements, computations can be
289: done while messages are in transition.
290: */
291: DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
292: DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
294: /*
295: Get pointers to vector data
296: */
297: DMDAVecGetArray(user->da,localX,&x);
298: DMDAVecGetArray(user->da,F,&f);
300: /*
301: Get local grid boundaries
302: */
303: DMDAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);
305: /*
306: Compute function over the locally owned part of the grid
307: */
308: for (k=zs; k<zs+zm; k++) {
309: for (j=ys; j<ys+ym; j++) {
310: for (i=xs; i<xs+xm; i++) {
311: if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
312: f[k][j][i] = x[k][j][i];
313: } else {
314: u = x[k][j][i];
315: u_east = x[k][j][i+1];
316: u_west = x[k][j][i-1];
317: u_north = x[k][j+1][i];
318: u_south = x[k][j-1][i];
319: u_up = x[k+1][j][i];
320: u_down = x[k-1][j][i];
321: u_xx = (-u_east + two*u - u_west)*hyhzdhx;
322: u_yy = (-u_north + two*u - u_south)*hxhzdhy;
323: u_zz = (-u_up + two*u - u_down)*hxhydhz;
324: f[k][j][i] = u_xx + u_yy + u_zz - sc*PetscExpScalar(u);
325: }
326: }
327: }
328: }
330: /*
331: Restore vectors
332: */
333: DMDAVecRestoreArray(user->da,localX,&x);
334: DMDAVecRestoreArray(user->da,F,&f);
335: DMRestoreLocalVector(user->da,&localX);
336: PetscLogFlops(11.0*ym*xm);
337: return(0);
338: }
339: /* ------------------------------------------------------------------- */
342: /*
343: FormJacobian - Evaluates Jacobian matrix.
345: Input Parameters:
346: . snes - the SNES context
347: . x - input vector
348: . ptr - optional user-defined context, as set by SNESSetJacobian()
350: Output Parameters:
351: . A - Jacobian matrix
352: . B - optionally different preconditioning matrix
353: . flag - flag indicating matrix structure
355: */
356: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
357: {
358: AppCtx *user = (AppCtx*)ptr; /* user-defined application context */
359: Mat jac = *B; /* Jacobian matrix */
360: Vec localX;
362: PetscInt i,j,k,Mx,My,Mz;
363: MatStencil col[7],row;
364: PetscInt xs,ys,zs,xm,ym,zm;
365: PetscScalar lambda,v[7],hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc,***x;
368: DMGetLocalVector(user->da,&localX);
369: DMDAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,
370: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
372: lambda = user->param;
373: hx = 1.0/(PetscReal)(Mx-1);
374: hy = 1.0/(PetscReal)(My-1);
375: hz = 1.0/(PetscReal)(Mz-1);
376: sc = hx*hy*hz*lambda;
377: hxhzdhy = hx*hz/hy;
378: hyhzdhx = hy*hz/hx;
379: hxhydhz = hx*hy/hz;
381: /*
382: Scatter ghost points to local vector, using the 2-step process
383: DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
384: By placing code between these two statements, computations can be
385: done while messages are in transition.
386: */
387: DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
388: DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
390: /*
391: Get pointer to vector data
392: */
393: DMDAVecGetArray(user->da,localX,&x);
395: /*
396: Get local grid boundaries
397: */
398: DMDAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);
400: /*
401: Compute entries for the locally owned part of the Jacobian.
402: - Currently, all PETSc parallel matrix formats are partitioned by
403: contiguous chunks of rows across the processors.
404: - Each processor needs to insert only elements that it owns
405: locally (but any non-local elements will be sent to the
406: appropriate processor during matrix assembly).
407: - Here, we set all entries for a particular row at once.
408: - We can set matrix entries either using either
409: MatSetValuesLocal() or MatSetValues(), as discussed above.
410: */
411: for (k=zs; k<zs+zm; k++) {
412: for (j=ys; j<ys+ym; j++) {
413: for (i=xs; i<xs+xm; i++) {
414: row.k = k; row.j = j; row.i = i;
415: /* boundary points */
416: if (i == 0 || j == 0 || k == 0|| i == Mx-1 || j == My-1 || k == Mz-1) {
417: v[0] = 1.0;
418: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
419: } else {
420: /* interior grid points */
421: v[0] = -hxhydhz; col[0].k=k-1;col[0].j=j; col[0].i = i;
422: v[1] = -hxhzdhy; col[1].k=k; col[1].j=j-1;col[1].i = i;
423: v[2] = -hyhzdhx; col[2].k=k; col[2].j=j; col[2].i = i-1;
424: 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;
425: v[4] = -hyhzdhx; col[4].k=k; col[4].j=j; col[4].i = i+1;
426: v[5] = -hxhzdhy; col[5].k=k; col[5].j=j+1;col[5].i = i;
427: v[6] = -hxhydhz; col[6].k=k+1;col[6].j=j; col[6].i = i;
428: MatSetValuesStencil(jac,1,&row,7,col,v,INSERT_VALUES);
429: }
430: }
431: }
432: }
433: DMDAVecRestoreArray(user->da,localX,&x);
434: DMRestoreLocalVector(user->da,&localX);
436: /*
437: Assemble matrix, using the 2-step process:
438: MatAssemblyBegin(), MatAssemblyEnd().
439: */
440: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
441: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
443: /*
444: Normally since the matrix has already been assembled above; this
445: would do nothing. But in the matrix free mode -snes_mf_operator
446: this tells the "matrix-free" matrix that a new linear system solve
447: is about to be done.
448: */
450: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
451: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
453: /*
454: Set flag to indicate that the Jacobian matrix retains an identical
455: nonzero structure throughout all nonlinear iterations (although the
456: values of the entries change). Thus, we can save some work in setting
457: up the preconditioner (e.g., no need to redo symbolic factorization for
458: ILU/ICC preconditioners).
459: - If the nonzero structure of the matrix is different during
460: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
461: must be used instead. If you are unsure whether the matrix
462: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
463: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
464: believes your assertion and does not check the structure
465: of the matrix. If you erroneously claim that the structure
466: is the same when it actually is not, the new preconditioner
467: will not function correctly. Thus, use this optimization
468: feature with caution!
469: */
470: *flag = SAME_NONZERO_PATTERN;
473: /*
474: Tell the matrix we will never add a new nonzero location to the
475: matrix. If we do, it will generate an error.
476: */
477: MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
478: return(0);
479: }