Actual source code: ex14.c
petsc-3.7.3 2016-08-01
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 <petscdm.h>
44: #include <petscdmda.h>
45: #include <petscsnes.h>
48: /*
49: User-defined application context - contains data needed by the
50: application-provided call-back routines, FormJacobian() and
51: FormFunction().
52: */
53: typedef struct {
54: PetscReal param; /* test problem parameter */
55: DM da; /* distributed array data structure */
56: } AppCtx;
58: /*
59: User-defined routines
60: */
61: extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*),FormInitialGuess(AppCtx*,Vec);
62: extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,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: PetscInt its; /* iterations for convergence */
73: MatFDColoring matfdcoloring;
74: PetscBool matrix_free = PETSC_FALSE,coloring = PETSC_FALSE, coloring_ds = PETSC_FALSE;
76: PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.,fnorm;
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(NULL,NULL,"-par",&user.param,NULL);
89: if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) SETERRQ(PETSC_COMM_SELF,1,"Lambda is out of range");
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Create nonlinear solver context
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: SNESCreate(PETSC_COMM_WORLD,&snes);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Create distributed array (DMDA) to manage parallel grid and vectors
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,-4,PETSC_DECIDE,PETSC_DECIDE,
100: PETSC_DECIDE,1,1,NULL,NULL,NULL,&user.da);
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Extract global vectors from DMDA; then duplicate for remaining
104: vectors that are the same types
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: DMCreateGlobalVector(user.da,&x);
107: VecDuplicate(x,&r);
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Set function evaluation routine and vector
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112: SNESSetFunction(snes,r,FormFunction,(void*)&user);
114: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: Create matrix data structure; set Jacobian evaluation routine
117: Set Jacobian matrix data structure and default Jacobian evaluation
118: routine. User can override with:
119: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
120: (unless user explicitly sets preconditioner)
121: -snes_mf_operator : form preconditioning matrix as set by the user,
122: but use matrix-free approx for Jacobian-vector
123: products within Newton-Krylov method
124: -fdcoloring : using finite differences with coloring to compute the Jacobian
126: Note one can use -matfd_coloring wp or ds the only reason for the -fdcoloring_ds option
127: below is to test the call to MatFDColoringSetType().
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129: PetscOptionsGetBool(NULL,NULL,"-snes_mf",&matrix_free,NULL);
130: PetscOptionsGetBool(NULL,NULL,"-fdcoloring",&coloring,NULL);
131: PetscOptionsGetBool(NULL,NULL,"-fdcoloring_ds",&coloring_ds,NULL);
132: if (!matrix_free) {
133: DMSetMatType(user.da,MATAIJ);
134: DMCreateMatrix(user.da,&J);
135: if (coloring) {
136: ISColoring iscoloring;
137: DMCreateColoring(user.da,IS_COLORING_GLOBAL,&iscoloring);
138: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
139: if (coloring_ds) {
140: MatFDColoringSetType(matfdcoloring,MATMFFD_DS);
141: }
142: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormFunction,&user);
143: MatFDColoringSetFromOptions(matfdcoloring);
144: MatFDColoringSetUp(J,iscoloring,matfdcoloring);
145: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
146: ISColoringDestroy(&iscoloring);
147: } else {
148: SNESSetJacobian(snes,J,J,FormJacobian,&user);
149: }
150: }
152: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153: Customize nonlinear solver; set runtime options
154: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155: SNESSetDM(snes,user.da);
156: SNESSetFromOptions(snes);
158: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159: Evaluate initial guess
160: Note: The user should initialize the vector, x, with the initial guess
161: for the nonlinear solver prior to calling SNESSolve(). In particular,
162: to employ an initial guess of zero, the user should explicitly set
163: this vector to zero by calling VecSet().
164: */
165: FormInitialGuess(&user,x);
167: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168: Solve nonlinear system
169: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170: SNESSolve(snes,NULL,x);
171: SNESGetIterationNumber(snes,&its);
173: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174: Explicitly check norm of the residual of the solution
175: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
176: FormFunction(snes,x,r,(void*)&user);
177: VecNorm(r,NORM_2,&fnorm);
178: PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D fnorm %g\n",its,(double)fnorm);
180: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181: Free work space. All PETSc objects should be destroyed when they
182: are no longer needed.
183: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185: if (!matrix_free) {
186: MatDestroy(&J);
187: }
188: VecDestroy(&x);
189: VecDestroy(&r);
190: SNESDestroy(&snes);
191: DMDestroy(&user.da);
192: if (coloring) {MatFDColoringDestroy(&matfdcoloring);}
193: PetscFinalize();
194: return(0);
195: }
196: /* ------------------------------------------------------------------- */
199: /*
200: FormInitialGuess - Forms initial approximation.
202: Input Parameters:
203: user - user-defined application context
204: X - vector
206: Output Parameter:
207: X - vector
208: */
209: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
210: {
211: PetscInt i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
213: PetscReal lambda,temp1,hx,hy,hz,tempk,tempj;
214: PetscScalar ***x;
217: 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);
219: lambda = user->param;
220: hx = 1.0/(PetscReal)(Mx-1);
221: hy = 1.0/(PetscReal)(My-1);
222: hz = 1.0/(PetscReal)(Mz-1);
223: temp1 = lambda/(lambda + 1.0);
225: /*
226: Get a pointer to vector data.
227: - For default PETSc vectors, VecGetArray() returns a pointer to
228: the data array. Otherwise, the routine is implementation dependent.
229: - You MUST call VecRestoreArray() when you no longer need access to
230: the array.
231: */
232: DMDAVecGetArray(user->da,X,&x);
234: /*
235: Get local grid boundaries (for 3-dimensional DMDA):
236: xs, ys, zs - starting grid indices (no ghost points)
237: xm, ym, zm - widths of local grid (no ghost points)
239: */
240: DMDAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);
242: /*
243: Compute initial guess over the locally owned part of the grid
244: */
245: for (k=zs; k<zs+zm; k++) {
246: tempk = (PetscReal)(PetscMin(k,Mz-k-1))*hz;
247: for (j=ys; j<ys+ym; j++) {
248: tempj = PetscMin((PetscReal)(PetscMin(j,My-j-1))*hy,tempk);
249: for (i=xs; i<xs+xm; i++) {
250: if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
251: /* boundary conditions are all zero Dirichlet */
252: x[k][j][i] = 0.0;
253: } else {
254: x[k][j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,tempj));
255: }
256: }
257: }
258: }
260: /*
261: Restore vector
262: */
263: DMDAVecRestoreArray(user->da,X,&x);
264: return(0);
265: }
266: /* ------------------------------------------------------------------- */
269: /*
270: FormFunction - Evaluates nonlinear function, F(x).
272: Input Parameters:
273: . snes - the SNES context
274: . X - input vector
275: . ptr - optional user-defined context, as set by SNESSetFunction()
277: Output Parameter:
278: . F - function vector
279: */
280: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
281: {
282: AppCtx *user = (AppCtx*)ptr;
284: PetscInt i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm;
285: PetscReal two = 2.0,lambda,hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc;
286: PetscScalar u_north,u_south,u_east,u_west,u_up,u_down,u,u_xx,u_yy,u_zz,***x,***f;
287: Vec localX;
288: DM da;
291: SNESGetDM(snes,&da);
292: DMGetLocalVector(da,&localX);
293: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
295: lambda = user->param;
296: hx = 1.0/(PetscReal)(Mx-1);
297: hy = 1.0/(PetscReal)(My-1);
298: hz = 1.0/(PetscReal)(Mz-1);
299: sc = hx*hy*hz*lambda;
300: hxhzdhy = hx*hz/hy;
301: hyhzdhx = hy*hz/hx;
302: hxhydhz = hx*hy/hz;
304: /*
305: Scatter ghost points to local vector,using the 2-step process
306: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
307: By placing code between these two statements, computations can be
308: done while messages are in transition.
309: */
310: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
311: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
313: /*
314: Get pointers to vector data
315: */
316: DMDAVecGetArrayRead(da,localX,&x);
317: DMDAVecGetArray(da,F,&f);
319: /*
320: Get local grid boundaries
321: */
322: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
324: /*
325: Compute function over the locally owned part of the grid
326: */
327: for (k=zs; k<zs+zm; k++) {
328: for (j=ys; j<ys+ym; j++) {
329: for (i=xs; i<xs+xm; i++) {
330: if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) {
331: f[k][j][i] = x[k][j][i];
332: } else {
333: u = x[k][j][i];
334: u_east = x[k][j][i+1];
335: u_west = x[k][j][i-1];
336: u_north = x[k][j+1][i];
337: u_south = x[k][j-1][i];
338: u_up = x[k+1][j][i];
339: u_down = x[k-1][j][i];
340: u_xx = (-u_east + two*u - u_west)*hyhzdhx;
341: u_yy = (-u_north + two*u - u_south)*hxhzdhy;
342: u_zz = (-u_up + two*u - u_down)*hxhydhz;
343: f[k][j][i] = u_xx + u_yy + u_zz - sc*PetscExpScalar(u);
344: }
345: }
346: }
347: }
349: /*
350: Restore vectors
351: */
352: DMDAVecRestoreArrayRead(da,localX,&x);
353: DMDAVecRestoreArray(da,F,&f);
354: DMRestoreLocalVector(da,&localX);
355: PetscLogFlops(11.0*ym*xm);
356: return(0);
357: }
358: /* ------------------------------------------------------------------- */
361: /*
362: FormJacobian - Evaluates Jacobian matrix.
364: Input Parameters:
365: . snes - the SNES context
366: . x - input vector
367: . ptr - optional user-defined context, as set by SNESSetJacobian()
369: Output Parameters:
370: . A - Jacobian matrix
371: . B - optionally different preconditioning matrix
373: */
374: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr)
375: {
376: AppCtx *user = (AppCtx*)ptr; /* user-defined application context */
377: Vec localX;
379: PetscInt i,j,k,Mx,My,Mz;
380: MatStencil col[7],row;
381: PetscInt xs,ys,zs,xm,ym,zm;
382: PetscScalar lambda,v[7],hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc,***x;
383: DM da;
386: SNESGetDM(snes,&da);
387: DMGetLocalVector(da,&localX);
388: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
390: lambda = user->param;
391: hx = 1.0/(PetscReal)(Mx-1);
392: hy = 1.0/(PetscReal)(My-1);
393: hz = 1.0/(PetscReal)(Mz-1);
394: sc = hx*hy*hz*lambda;
395: hxhzdhy = hx*hz/hy;
396: hyhzdhx = hy*hz/hx;
397: hxhydhz = hx*hy/hz;
399: /*
400: Scatter ghost points to local vector, using the 2-step process
401: DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
402: By placing code between these two statements, computations can be
403: done while messages are in transition.
404: */
405: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
406: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
408: /*
409: Get pointer to vector data
410: */
411: DMDAVecGetArrayRead(da,localX,&x);
413: /*
414: Get local grid boundaries
415: */
416: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
418: /*
419: Compute entries for the locally owned part of the Jacobian.
420: - Currently, all PETSc parallel matrix formats are partitioned by
421: contiguous chunks of rows across the processors.
422: - Each processor needs to insert only elements that it owns
423: locally (but any non-local elements will be sent to the
424: appropriate processor during matrix assembly).
425: - Here, we set all entries for a particular row at once.
426: - We can set matrix entries either using either
427: MatSetValuesLocal() or MatSetValues(), as discussed above.
428: */
429: for (k=zs; k<zs+zm; k++) {
430: for (j=ys; j<ys+ym; j++) {
431: for (i=xs; i<xs+xm; i++) {
432: row.k = k; row.j = j; row.i = i;
433: /* boundary points */
434: if (i == 0 || j == 0 || k == 0|| i == Mx-1 || j == My-1 || k == Mz-1) {
435: v[0] = 1.0;
436: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
437: } else {
438: /* interior grid points */
439: v[0] = -hxhydhz; col[0].k=k-1;col[0].j=j; col[0].i = i;
440: v[1] = -hxhzdhy; col[1].k=k; col[1].j=j-1;col[1].i = i;
441: v[2] = -hyhzdhx; col[2].k=k; col[2].j=j; col[2].i = i-1;
442: 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;
443: v[4] = -hyhzdhx; col[4].k=k; col[4].j=j; col[4].i = i+1;
444: v[5] = -hxhzdhy; col[5].k=k; col[5].j=j+1;col[5].i = i;
445: v[6] = -hxhydhz; col[6].k=k+1;col[6].j=j; col[6].i = i;
446: MatSetValuesStencil(jac,1,&row,7,col,v,INSERT_VALUES);
447: }
448: }
449: }
450: }
451: DMDAVecRestoreArrayRead(da,localX,&x);
452: DMRestoreLocalVector(da,&localX);
454: /*
455: Assemble matrix, using the 2-step process:
456: MatAssemblyBegin(), MatAssemblyEnd().
457: */
458: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
459: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
461: /*
462: Normally since the matrix has already been assembled above; this
463: would do nothing. But in the matrix free mode -snes_mf_operator
464: this tells the "matrix-free" matrix that a new linear system solve
465: is about to be done.
466: */
468: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
469: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
471: /*
472: Tell the matrix we will never add a new nonzero location to the
473: matrix. If we do, it will generate an error.
474: */
475: MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
476: return(0);
477: }