Actual source code: ex5.c
2: /* Program usage: mpirun -np <procs> ex5 [-help] [all PETSc options] */
4: static char help[] = "Bratu nonlinear PDE in 2d.\n\
5: We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\
6: domain, using distributed arrays (DAs) to partition the parallel grid.\n\
7: The command line options include:\n\
8: -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
9: problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n";
11: /*T
12: Concepts: SNES^parallel Bratu example
13: Concepts: DA^using distributed arrays;
14: Processors: n
15: T*/
17: /* ------------------------------------------------------------------------
19: Solid Fuel Ignition (SFI) problem. This problem is modeled by
20: the partial differential equation
21:
22: -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
23:
24: with boundary conditions
25:
26: u = 0 for x = 0, x = 1, y = 0, y = 1.
27:
28: A finite difference approximation with the usual 5-point stencil
29: is used to discretize the boundary value problem to obtain a nonlinear
30: system of equations.
33: ------------------------------------------------------------------------- */
35: /*
36: Include "petscda.h" so that we can use distributed arrays (DAs).
37: Include "petscsnes.h" so that we can use SNES solvers. Note that this
38: file automatically includes:
39: petsc.h - base PETSc routines petscvec.h - vectors
40: petscsys.h - system routines petscmat.h - matrices
41: petscis.h - index sets petscksp.h - Krylov subspace methods
42: petscviewer.h - viewers petscpc.h - preconditioners
43: petscksp.h - linear solvers
44: */
45: #include petscda.h
46: #include petscsnes.h
48: /*
49: User-defined application context - contains data needed by the
50: application-provided call-back routines, FormJacobianLocal() and
51: FormFunctionLocal().
52: */
53: typedef struct {
54: DA da; /* distributed array data structure */
55: PassiveReal param; /* test problem parameter */
56: } AppCtx;
58: /*
59: User-defined routines
60: */
67: int main(int argc,char **argv)
68: {
69: SNES snes; /* nonlinear solver */
70: Vec x,r; /* solution, residual vectors */
71: Mat A,J; /* Jacobian matrix */
72: AppCtx user; /* user-defined work context */
73: PetscInt its; /* iterations for convergence */
74: PetscTruth matlab_function = PETSC_FALSE;
75: PetscTruth fd_jacobian = PETSC_FALSE,adic_jacobian=PETSC_FALSE;
76: PetscTruth adicmf_jacobian = PETSC_FALSE;
77: PetscErrorCode ierr;
78: PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
79: MatFDColoring matfdcoloring = 0;
80: ISColoring iscoloring;
83:
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Initialize program
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: PetscInitialize(&argc,&argv,(char *)0,help);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Initialize problem parameters
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: user.param = 6.0;
94: PetscOptionsGetReal(PETSC_NULL,"-par",&user.param,PETSC_NULL);
95: if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
96: SETERRQ(1,"Lambda is out of range");
97: }
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Create nonlinear solver context
101: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102: SNESCreate(PETSC_COMM_WORLD,&snes);
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Create distributed array (DA) to manage parallel grid and vectors
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,
108: 1,1,PETSC_NULL,PETSC_NULL,&user.da);
110: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: Extract global vectors from DA; then duplicate for remaining
112: vectors that are the same types
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: DACreateGlobalVector(user.da,&x);
115: VecDuplicate(x,&r);
117: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118: Set local function evaluation routine
119: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120: DASetLocalFunction(user.da,(DALocalFunction1)FormFunctionLocal);
121: DASetLocalJacobian(user.da,(DALocalFunction1)FormJacobianLocal);
122: DASetLocalAdicFunction(user.da,ad_FormFunctionLocal);
124: /* Decide which FormFunction to use */
125: PetscOptionsGetLogical(PETSC_NULL,"-matlab_function",&matlab_function,0);
127: SNESSetFunction(snes,r,SNESDAFormFunction,&user);
128: #if defined(PETSC_HAVE_MATLAB) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
129: if (matlab_function) {
130: SNESSetFunction(snes,r,FormFunctionMatlab,&user);
131: }
132: #endif
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Create matrix data structure; set Jacobian evaluation routine
137: Set Jacobian matrix data structure and default Jacobian evaluation
138: routine. User can override with:
139: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
140: (unless user explicitly sets preconditioner)
141: -snes_mf_operator : form preconditioning matrix as set by the user,
142: but use matrix-free approx for Jacobian-vector
143: products within Newton-Krylov method
145: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146: DAGetMatrix(user.da,MATAIJ,&J);
147: A = J;
149: PetscOptionsGetLogical(PETSC_NULL,"-fd_jacobian",&fd_jacobian,0);
150: PetscOptionsGetLogical(PETSC_NULL,"-adic_jacobian",&adic_jacobian,0);
151: PetscOptionsGetLogical(PETSC_NULL,"-adicmf_jacobian",&adicmf_jacobian,0);
152: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
153: if (adicmf_jacobian) {
154: DASetLocalAdicMFFunction(user.da,admf_FormFunctionLocal);
155: MatRegisterDAAD();
156: MatCreateDAAD(user.da,&A);
157: MatDAADSetSNES(A,snes);
158: MatDAADSetCtx(A,&user);
159: }
160: #endif
162: if (fd_jacobian) {
163: DAGetColoring(user.da,IS_COLORING_LOCAL,&iscoloring);
164: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
165: ISColoringDestroy(iscoloring);
166: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESDAFormFunction,&user);
167: MatFDColoringSetFromOptions(matfdcoloring);
168: SNESSetJacobian(snes,A,J,SNESDefaultComputeJacobianColor,matfdcoloring);
169: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
170: } else if (adic_jacobian) {
171: DAGetColoring(user.da,IS_COLORING_GHOSTED,&iscoloring);
172: MatSetColoring(J,iscoloring);
173: ISColoringDestroy(iscoloring);
174: SNESSetJacobian(snes,A,J,SNESDAComputeJacobianWithAdic,&user);
175: #endif
176: } else {
177: SNESSetJacobian(snes,A,J,SNESDAComputeJacobian,&user);
178: }
180: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181: Customize nonlinear solver; set runtime options
182: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
183: SNESSetFromOptions(snes);
185: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186: Evaluate initial guess
187: Note: The user should initialize the vector, x, with the initial guess
188: for the nonlinear solver prior to calling SNESSolve(). In particular,
189: to employ an initial guess of zero, the user should explicitly set
190: this vector to zero by calling VecSet().
191: */
192: FormInitialGuess(&user,x);
194: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195: Solve nonlinear system
196: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197: SNESSolve(snes,x);
198: SNESGetIterationNumber(snes,&its);
200: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %D\n",its);
204: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205: Free work space. All PETSc objects should be destroyed when they
206: are no longer needed.
207: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209: if (A != J) {
210: MatDestroy(A);
211: }
212: MatDestroy(J);
213: if (matfdcoloring) {
214: MatFDColoringDestroy(matfdcoloring);
215: }
216: VecDestroy(x);
217: VecDestroy(r);
218: SNESDestroy(snes);
219: DADestroy(user.da);
220: PetscFinalize();
222: return(0);
223: }
224: /* ------------------------------------------------------------------- */
227: /*
228: FormInitialGuess - Forms initial approximation.
230: Input Parameters:
231: user - user-defined application context
232: X - vector
234: Output Parameter:
235: X - vector
236: */
237: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
238: {
239: PetscInt i,j,Mx,My,xs,ys,xm,ym;
241: PetscReal lambda,temp1,temp,hx,hy;
242: PetscScalar **x;
245: DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
246: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
248: lambda = user->param;
249: hx = 1.0/(PetscReal)(Mx-1);
250: hy = 1.0/(PetscReal)(My-1);
251: temp1 = lambda/(lambda + 1.0);
253: /*
254: Get a pointer to vector data.
255: - For default PETSc vectors, VecGetArray() returns a pointer to
256: the data array. Otherwise, the routine is implementation dependent.
257: - You MUST call VecRestoreArray() when you no longer need access to
258: the array.
259: */
260: DAVecGetArray(user->da,X,&x);
262: /*
263: Get local grid boundaries (for 2-dimensional DA):
264: xs, ys - starting grid indices (no ghost points)
265: xm, ym - widths of local grid (no ghost points)
267: */
268: DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
270: /*
271: Compute initial guess over the locally owned part of the grid
272: */
273: for (j=ys; j<ys+ym; j++) {
274: temp = (PetscReal)(PetscMin(j,My-j-1))*hy;
275: for (i=xs; i<xs+xm; i++) {
277: if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
278: /* boundary conditions are all zero Dirichlet */
279: x[j][i] = 0.0;
280: } else {
281: x[j][i] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
282: }
283: }
284: }
286: /*
287: Restore vector
288: */
289: DAVecRestoreArray(user->da,X,&x);
291: return(0);
292: }
293: /* ------------------------------------------------------------------- */
296: /*
297: FormFunctionLocal - Evaluates nonlinear function, F(x).
299: Process adiC(36): FormFunctionLocal
301: */
302: PetscErrorCode FormFunctionLocal(DALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
303: {
305: PetscInt i,j;
306: PetscReal two = 2.0,lambda,hx,hy,hxdhy,hydhx,sc;
307: PetscScalar u,uxx,uyy;
311: lambda = user->param;
312: hx = 1.0/(PetscReal)(info->mx-1);
313: hy = 1.0/(PetscReal)(info->my-1);
314: sc = hx*hy*lambda;
315: hxdhy = hx/hy;
316: hydhx = hy/hx;
317: /*
318: Compute function over the locally owned part of the grid
319: */
320: for (j=info->ys; j<info->ys+info->ym; j++) {
321: for (i=info->xs; i<info->xs+info->xm; i++) {
322: if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
323: f[j][i] = x[j][i];
324: } else {
325: u = x[j][i];
326: uxx = (two*u - x[j][i-1] - x[j][i+1])*hydhx;
327: uyy = (two*u - x[j-1][i] - x[j+1][i])*hxdhy;
328: f[j][i] = uxx + uyy - sc*PetscExpScalar(u);
329: }
330: }
331: }
333: PetscLogFlops(11*info->ym*info->xm);
334: return(0);
335: }
339: /*
340: FormJacobianLocal - Evaluates Jacobian matrix.
343: */
344: PetscErrorCode FormJacobianLocal(DALocalInfo *info,PetscScalar **x,Mat jac,AppCtx *user)
345: {
347: PetscInt i,j;
348: MatStencil col[5],row;
349: PetscScalar lambda,v[5],hx,hy,hxdhy,hydhx,sc;
352: lambda = user->param;
353: hx = 1.0/(PetscReal)(info->mx-1);
354: hy = 1.0/(PetscReal)(info->my-1);
355: sc = hx*hy*lambda;
356: hxdhy = hx/hy;
357: hydhx = hy/hx;
360: /*
361: Compute entries for the locally owned part of the Jacobian.
362: - Currently, all PETSc parallel matrix formats are partitioned by
363: contiguous chunks of rows across the processors.
364: - Each processor needs to insert only elements that it owns
365: locally (but any non-local elements will be sent to the
366: appropriate processor during matrix assembly).
367: - Here, we set all entries for a particular row at once.
368: - We can set matrix entries either using either
369: MatSetValuesLocal() or MatSetValues(), as discussed above.
370: */
371: for (j=info->ys; j<info->ys+info->ym; j++) {
372: for (i=info->xs; i<info->xs+info->xm; i++) {
373: row.j = j; row.i = i;
374: /* boundary points */
375: if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
376: v[0] = 1.0;
377: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
378: } else {
379: /* interior grid points */
380: v[0] = -hxdhy; col[0].j = j - 1; col[0].i = i;
381: v[1] = -hydhx; col[1].j = j; col[1].i = i-1;
382: v[2] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(x[j][i]); col[2].j = row.j; col[2].i = row.i;
383: v[3] = -hydhx; col[3].j = j; col[3].i = i+1;
384: v[4] = -hxdhy; col[4].j = j + 1; col[4].i = i;
385: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
386: }
387: }
388: }
390: /*
391: Assemble matrix, using the 2-step process:
392: MatAssemblyBegin(), MatAssemblyEnd().
393: */
394: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
395: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
396: /*
397: Tell the matrix we will never add a new nonzero location to the
398: matrix. If we do, it will generate an error.
399: */
400: MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR);
401: return(0);
402: }
404: /*
405: Variant of FormFunction() that computes the function in Matlab
406: */
407: #if defined(PETSC_HAVE_MATLAB) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
408: PetscErrorCode FormFunctionMatlab(SNES snes,Vec X,Vec F,void *ptr)
409: {
410: AppCtx *user = (AppCtx*)ptr;
412: PetscInt Mx,My;
413: PetscReal lambda,hx,hy;
414: Vec localX,localF;
415: MPI_Comm comm;
418: DAGetLocalVector(user->da,&localX);
419: DAGetLocalVector(user->da,&localF);
420: PetscObjectSetName((PetscObject)localX,"localX");
421: PetscObjectSetName((PetscObject)localF,"localF");
422: DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
423: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
425: lambda = user->param;
426: hx = 1.0/(PetscReal)(Mx-1);
427: hy = 1.0/(PetscReal)(My-1);
429: PetscObjectGetComm((PetscObject)snes,&comm);
430: /*
431: Scatter ghost points to local vector,using the 2-step process
432: DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
433: By placing code between these two statements, computations can be
434: done while messages are in transition.
435: */
436: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
437: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
438: PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localX);
439: PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm),"localF=ex5m(localX,%18.16e,%18.16e,%18.16e)",hx,hy,lambda);
440: PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localF);
442: /*
443: Insert values into global vector
444: */
445: DALocalToGlobal(user->da,localF,INSERT_VALUES,F);
446: DARestoreLocalVector(user->da,&localX);
447: DARestoreLocalVector(user->da,&localF);
448: return(0);
449: }
450: #endif