Actual source code: ex5.c
2: static char help[] = "Bratu nonlinear PDE in 2d.\n\
3: We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\
4: domain, using distributed arrays (DAs) 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: DA^using distributed arrays;
12: Concepts: IS coloirng types;
13: Processors: n
14: T*/
16: /* ------------------------------------------------------------------------
18: Solid Fuel Ignition (SFI) problem. This problem is modeled by
19: the partial differential equation
20:
21: -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
22:
23: with boundary conditions
24:
25: u = 0 for x = 0, x = 1, y = 0, y = 1.
26:
27: A finite difference approximation with the usual 5-point stencil
28: is used to discretize the boundary value problem to obtain a nonlinear
29: system of equations.
31: Program usage: mpiexec -n <procs> ex5 [-help] [all PETSc options]
32: e.g.,
33: ./ex5 -fd_jacobian -mat_fd_coloring_view_draw -draw_pause -1
34: mpiexec -n 2 ./ex5 -fd_jacobian_ghosted -log_summary
36: ------------------------------------------------------------------------- */
38: /*
39: Include "petscda.h" so that we can use distributed arrays (DAs).
40: Include "petscsnes.h" so that we can use SNES solvers. Note that this
41: file automatically includes:
42: petsc.h - base PETSc routines petscvec.h - vectors
43: petscsys.h - system routines petscmat.h - matrices
44: petscis.h - index sets petscksp.h - Krylov subspace methods
45: petscviewer.h - viewers petscpc.h - preconditioners
46: petscksp.h - linear solvers
47: */
48: #include petscda.h
49: #include petscsnes.h
51: /*
52: User-defined application context - contains data needed by the
53: application-provided call-back routines, FormJacobianLocal() and
54: FormFunctionLocal().
55: */
56: typedef struct {
57: DA da; /* distributed array data structure */
58: PassiveReal param; /* test problem parameter */
59: } AppCtx;
61: /*
62: User-defined routines
63: */
71: int main(int argc,char **argv)
72: {
73: SNES snes; /* nonlinear solver */
74: Vec x,r; /* solution, residual vectors */
75: Mat A,J; /* Jacobian matrix */
76: AppCtx user; /* user-defined work context */
77: PetscInt its; /* iterations for convergence */
78: PetscTruth matlab_function = PETSC_FALSE;
79: PetscTruth fd_jacobian = PETSC_FALSE,adic_jacobian=PETSC_FALSE,fd_jacobian_ghosted=PETSC_FALSE;
80: PetscTruth adicmf_jacobian = PETSC_FALSE;
81: PetscErrorCode ierr;
82: PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
83: MatFDColoring matfdcoloring = 0;
84: ISColoring iscoloring;
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Initialize program
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: PetscInitialize(&argc,&argv,(char *)0,help);
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Initialize problem parameters
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: user.param = 6.0;
96: PetscOptionsGetReal(PETSC_NULL,"-par",&user.param,PETSC_NULL);
97: if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
98: SETERRQ(1,"Lambda is out of range");
99: }
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Create nonlinear solver context
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: SNESCreate(PETSC_COMM_WORLD,&snes);
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Create distributed array (DA) to manage parallel grid and vectors
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,
110: 1,1,PETSC_NULL,PETSC_NULL,&user.da);
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Extract global vectors from DA; then duplicate for remaining
114: vectors that are the same types
115: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: DACreateGlobalVector(user.da,&x);
117: VecDuplicate(x,&r);
119: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: Create matrix data structure; set Jacobian evaluation routine
122: Set Jacobian matrix data structure and default Jacobian evaluation
123: routine. User can override with:
124: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
125: (unless user explicitly sets preconditioner)
126: -snes_mf_operator : form preconditioning matrix as set by the user,
127: but use matrix-free approx for Jacobian-vector
128: products within Newton-Krylov method
130: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131: /* J can be type of MATAIJ,MATBAIJ or MATSBAIJ */
132: DAGetMatrix(user.da,MATAIJ,&J);
133:
134: A = J;
135: PetscOptionsGetTruth(PETSC_NULL,"-fd_jacobian",&fd_jacobian,0);
136: PetscOptionsGetTruth(PETSC_NULL,"-fd_jacobian_ghosted",&fd_jacobian_ghosted,0);
137: PetscOptionsGetTruth(PETSC_NULL,"-adic_jacobian",&adic_jacobian,0);
138: PetscOptionsGetTruth(PETSC_NULL,"-adicmf_jacobian",&adicmf_jacobian,0);
139: #if defined(PETSC_HAVE_ADIC)
140: if (adicmf_jacobian) {
141: DASetLocalAdicMFFunction(user.da,admf_FormFunctionLocal);
142: MatRegisterDAAD();
143: MatCreateDAAD(user.da,&A);
144: MatDAADSetSNES(A,snes);
145: MatDAADSetCtx(A,&user);
146: }
147: #endif
149: if (fd_jacobian) {
150: DAGetColoring(user.da,IS_COLORING_GLOBAL,&iscoloring);
151: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
152: ISColoringDestroy(iscoloring);
153: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESDAFormFunction,&user);
154: MatFDColoringSetFromOptions(matfdcoloring);
155: SNESSetJacobian(snes,A,J,SNESDefaultComputeJacobianColor,matfdcoloring);
156: } else if (fd_jacobian_ghosted) {
157: DAGetColoring(user.da,IS_COLORING_GHOSTED,&iscoloring);
158: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
159: ISColoringDestroy(iscoloring);
160: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESDAFormFunction,&user);
161: MatFDColoringSetFromOptions(matfdcoloring);
162: /* now, call a customized SNESDefaultComputeJacobianColor() */
163: SNESSetJacobian(snes,A,J,MySNESDefaultComputeJacobianColor,matfdcoloring);
164: #if defined(PETSC_HAVE_ADIC)
165: } else if (adic_jacobian) {
166: DAGetColoring(user.da,IS_COLORING_GHOSTED,&iscoloring);
167: MatSetColoring(J,iscoloring);
168: ISColoringDestroy(iscoloring);
169: SNESSetJacobian(snes,A,J,SNESDAComputeJacobianWithAdic,&user);
170: #endif
171: } else {
172: SNESSetJacobian(snes,A,J,SNESDAComputeJacobian,&user);
173: }
175: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176: Set local function evaluation routine
177: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178: DASetLocalFunction(user.da,(DALocalFunction1)FormFunctionLocal);
179: DASetLocalJacobian(user.da,(DALocalFunction1)FormJacobianLocal);
180: DASetLocalAdicFunction(user.da,ad_FormFunctionLocal);
182: /* Decide which FormFunction to use */
183: PetscOptionsGetTruth(PETSC_NULL,"-matlab_function",&matlab_function,0);
185: SNESSetFunction(snes,r,SNESDAFormFunction,&user);
186: #if defined(PETSC_HAVE_MATLAB_ENGINE)
187: if (matlab_function) {
188: SNESSetFunction(snes,r,FormFunctionMatlab,&user);
189: }
190: #endif
192: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193: Customize nonlinear solver; set runtime options
194: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195: SNESSetFromOptions(snes);
197: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198: Evaluate initial guess
199: Note: The user should initialize the vector, x, with the initial guess
200: for the nonlinear solver prior to calling SNESSolve(). In particular,
201: to employ an initial guess of zero, the user should explicitly set
202: this vector to zero by calling VecSet().
203: */
205: {
206: PetscTruth test_appctx = PETSC_FALSE;
207: PetscOptionsGetTruth(PETSC_NULL,"-test_appctx",&test_appctx,0);
208: if (test_appctx) {
209: AppCtx *puser;
210: SNESSetApplicationContext(snes,&user);
211: SNESGetApplicationContext(snes,(void **)&puser);
212: FormInitialGuess(puser,x);
213: } else {
214: FormInitialGuess(&user,x);
215: }
216: }
218: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219: Solve nonlinear system
220: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
221: SNESSolve(snes,PETSC_NULL,x);
222: SNESGetIterationNumber(snes,&its);
224: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
226: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %D\n",its);
228: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229: Free work space. All PETSc objects should be destroyed when they
230: are no longer needed.
231: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
233: if (A != J) {
234: MatDestroy(A);
235: }
236: MatDestroy(J);
237: if (matfdcoloring) {
238: MatFDColoringDestroy(matfdcoloring);
239: }
240: VecDestroy(x);
241: VecDestroy(r);
242: SNESDestroy(snes);
243: DADestroy(user.da);
244: PetscFinalize();
246: return(0);
247: }
248: /* ------------------------------------------------------------------- */
251: /*
252: FormInitialGuess - Forms initial approximation.
254: Input Parameters:
255: user - user-defined application context
256: X - vector
258: Output Parameter:
259: X - vector
260: */
261: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
262: {
263: PetscInt i,j,Mx,My,xs,ys,xm,ym;
265: PetscReal lambda,temp1,temp,hx,hy;
266: PetscScalar **x;
269: DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
270: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
272: lambda = user->param;
273: hx = 1.0/(PetscReal)(Mx-1);
274: hy = 1.0/(PetscReal)(My-1);
275: temp1 = lambda/(lambda + 1.0);
277: /*
278: Get a pointer to vector data.
279: - For default PETSc vectors, VecGetArray() returns a pointer to
280: the data array. Otherwise, the routine is implementation dependent.
281: - You MUST call VecRestoreArray() when you no longer need access to
282: the array.
283: */
284: DAVecGetArray(user->da,X,&x);
286: /*
287: Get local grid boundaries (for 2-dimensional DA):
288: xs, ys - starting grid indices (no ghost points)
289: xm, ym - widths of local grid (no ghost points)
291: */
292: DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
294: /*
295: Compute initial guess over the locally owned part of the grid
296: */
297: for (j=ys; j<ys+ym; j++) {
298: temp = (PetscReal)(PetscMin(j,My-j-1))*hy;
299: for (i=xs; i<xs+xm; i++) {
300: if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
301: /* boundary conditions are all zero Dirichlet */
302: x[j][i] = 0.0;
303: } else {
304: x[j][i] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
305: }
306: }
307: }
309: /*
310: Restore vector
311: */
312: DAVecRestoreArray(user->da,X,&x);
314: return(0);
315: }
316: /* ------------------------------------------------------------------- */
319: /*
320: FormFunctionLocal - Evaluates nonlinear function, F(x).
322: Process adiC(36): FormFunctionLocal
324: */
325: PetscErrorCode FormFunctionLocal(DALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
326: {
328: PetscInt i,j;
329: PetscReal lambda,hx,hy,hxdhy,hydhx,sc;
330: PetscScalar u,uxx,uyy;
334: lambda = user->param;
335: hx = 1.0/(PetscReal)(info->mx-1);
336: hy = 1.0/(PetscReal)(info->my-1);
337: sc = hx*hy*lambda;
338: hxdhy = hx/hy;
339: hydhx = hy/hx;
340: /*
341: Compute function over the locally owned part of the grid
342: */
343: for (j=info->ys; j<info->ys+info->ym; j++) {
344: for (i=info->xs; i<info->xs+info->xm; i++) {
345: if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
346: f[j][i] = x[j][i];
347: } else {
348: u = x[j][i];
349: uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx;
350: uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy;
351: f[j][i] = uxx + uyy - sc*PetscExpScalar(u);
352: }
353: }
354: }
356: PetscLogFlops(11*info->ym*info->xm);
357: return(0);
358: }
362: /*
363: FormJacobianLocal - Evaluates Jacobian matrix.
364: */
365: PetscErrorCode FormJacobianLocal(DALocalInfo *info,PetscScalar **x,Mat jac,AppCtx *user)
366: {
368: PetscInt i,j;
369: MatStencil col[5],row;
370: PetscScalar lambda,v[5],hx,hy,hxdhy,hydhx,sc;
373: lambda = user->param;
374: hx = 1.0/(PetscReal)(info->mx-1);
375: hy = 1.0/(PetscReal)(info->my-1);
376: sc = hx*hy*lambda;
377: hxdhy = hx/hy;
378: hydhx = hy/hx;
381: /*
382: Compute entries for the locally owned part of the Jacobian.
383: - Currently, all PETSc parallel matrix formats are partitioned by
384: contiguous chunks of rows across the processors.
385: - Each processor needs to insert only elements that it owns
386: locally (but any non-local elements will be sent to the
387: appropriate processor during matrix assembly).
388: - Here, we set all entries for a particular row at once.
389: - We can set matrix entries either using either
390: MatSetValuesLocal() or MatSetValues(), as discussed above.
391: */
392: for (j=info->ys; j<info->ys+info->ym; j++) {
393: for (i=info->xs; i<info->xs+info->xm; i++) {
394: row.j = j; row.i = i;
395: /* boundary points */
396: if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
397: v[0] = 1.0;
398: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
399: } else {
400: /* interior grid points */
401: v[0] = -hxdhy; col[0].j = j - 1; col[0].i = i;
402: v[1] = -hydhx; col[1].j = j; col[1].i = i-1;
403: v[2] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(x[j][i]); col[2].j = row.j; col[2].i = row.i;
404: v[3] = -hydhx; col[3].j = j; col[3].i = i+1;
405: v[4] = -hxdhy; col[4].j = j + 1; col[4].i = i;
406: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
407: }
408: }
409: }
411: /*
412: Assemble matrix, using the 2-step process:
413: MatAssemblyBegin(), MatAssemblyEnd().
414: */
415: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
416: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
417: /*
418: Tell the matrix we will never add a new nonzero location to the
419: matrix. If we do, it will generate an error.
420: */
421: MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
422: return(0);
423: }
425: /*
426: Variant of FormFunction() that computes the function in Matlab
427: */
428: #if defined(PETSC_HAVE_MATLAB_ENGINE)
429: PetscErrorCode FormFunctionMatlab(SNES snes,Vec X,Vec F,void *ptr)
430: {
431: AppCtx *user = (AppCtx*)ptr;
433: PetscInt Mx,My;
434: PetscReal lambda,hx,hy;
435: Vec localX,localF;
436: MPI_Comm comm;
439: DAGetLocalVector(user->da,&localX);
440: DAGetLocalVector(user->da,&localF);
441: PetscObjectSetName((PetscObject)localX,"localX");
442: PetscObjectSetName((PetscObject)localF,"localF");
443: DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
444: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
446: lambda = user->param;
447: hx = 1.0/(PetscReal)(Mx-1);
448: hy = 1.0/(PetscReal)(My-1);
450: PetscObjectGetComm((PetscObject)snes,&comm);
451: /*
452: Scatter ghost points to local vector,using the 2-step process
453: DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
454: By placing code between these two statements, computations can be
455: done while messages are in transition.
456: */
457: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
458: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
459: PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localX);
460: PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm),"localF=ex5m(localX,%18.16e,%18.16e,%18.16e)",hx,hy,lambda);
461: PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localF);
463: /*
464: Insert values into global vector
465: */
466: DALocalToGlobal(user->da,localF,INSERT_VALUES,F);
467: DARestoreLocalVector(user->da,&localX);
468: DARestoreLocalVector(user->da,&localF);
469: return(0);
470: }
471: #endif
475: /*
476: MySNESDefaultComputeJacobianColor - Computes the Jacobian using
477: finite differences and coloring to exploit matrix sparsity.
478: It is customized from SNESDefaultComputeJacobianColor.
479: The input global vector x1 is scattered to x1_local
480: which then is passed into MatFDColoringApply() for reducing the
481: VecScatterBingin/End.
482: */
483: PetscErrorCode MySNESDefaultComputeJacobianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
484: {
485: MatFDColoring color = (MatFDColoring) ctx;
487: Vec f;
488: PetscErrorCode (*ff)(void),(*fd)(void);
489: void *fctx;
490: DA da;
491: Vec x1_loc;
494: *flag = SAME_NONZERO_PATTERN;
495: SNESGetFunction(snes,&f,(PetscErrorCode (**)(SNES,Vec,Vec,void*))&ff,0);
496: MatFDColoringGetFunction(color,&fd,&fctx);
497: if (fd == ff) { /* reuse function value computed in SNES */
498: MatFDColoringSetF(color,f);
499: }
500: /* Now, get x1_loc and scatter global x1 onto x1_loc */
501: da = *(DA*)fctx;
502: DAGetLocalVector(da,&x1_loc);
503: DAGlobalToLocalBegin(da,x1,INSERT_VALUES,x1_loc);
504: DAGlobalToLocalEnd(da,x1,INSERT_VALUES,x1_loc);
505: MatFDColoringApply(*B,color,x1_loc,flag,snes);
506: DARestoreLocalVector(da,&x1_loc);
507: if (*J != *B) {
508: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
509: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
510: }
511: return(0);
512: }