Actual source code: ex14.c
2: static char help[] = "Solves a nonlinear system in parallel with a user-defined Newton method.\n\
3: Uses KSP to solve the linearized Newton systems. This solver\n\
4: is a very simplistic inexact Newton method. The intent of this code is to\n\
5: demonstrate the repeated solution of linear systems with the same nonzero pattern.\n\
6: \n\
7: This is NOT the recommended approach for solving nonlinear problems with PETSc!\n\
8: We urge users to employ the SNES component for solving nonlinear problems whenever\n\
9: possible, as it offers many advantages over coding nonlinear solvers independently.\n\
10: \n\
11: We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\
12: domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\
13: The command line options include:\n\
14: -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
15: problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\
16: -mx <xg>, where <xg> = number of grid points in the x-direction\n\
17: -my <yg>, where <yg> = number of grid points in the y-direction\n\
18: -Nx <npx>, where <npx> = number of processors in the x-direction\n\
19: -Ny <npy>, where <npy> = number of processors in the y-direction\n\n";
21: /*T
22: Concepts: KSP^writing a user-defined nonlinear solver (parallel Bratu example);
23: Concepts: DMDA^using distributed arrays;
24: Processors: n
25: T*/
27: /* ------------------------------------------------------------------------
29: Solid Fuel Ignition (SFI) problem. This problem is modeled by
30: the partial differential equation
32: -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
34: with boundary conditions
36: u = 0 for x = 0, x = 1, y = 0, y = 1.
38: A finite difference approximation with the usual 5-point stencil
39: is used to discretize the boundary value problem to obtain a nonlinear
40: system of equations.
42: The SNES version of this problem is: snes/tutorials/ex5.c
43: We urge users to employ the SNES component for solving nonlinear
44: problems whenever possible, as it offers many advantages over coding
45: nonlinear solvers independently.
47: ------------------------------------------------------------------------- */
49: /*
50: Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
51: Include "petscksp.h" so that we can use KSP solvers. Note that this
52: file automatically includes:
53: petscsys.h - base PETSc routines petscvec.h - vectors
54: petscmat.h - matrices
55: petscis.h - index sets petscksp.h - Krylov subspace methods
56: petscviewer.h - viewers petscpc.h - preconditioners
57: */
58: #include <petscdm.h>
59: #include <petscdmda.h>
60: #include <petscksp.h>
62: /*
63: User-defined application context - contains data needed by the
64: application-provided call-back routines, ComputeJacobian() and
65: ComputeFunction().
66: */
67: typedef struct {
68: PetscReal param; /* test problem parameter */
69: PetscInt mx,my; /* discretization in x,y directions */
70: Vec localX; /* ghosted local vector */
71: DM da; /* distributed array data structure */
72: } AppCtx;
74: /*
75: User-defined routines
76: */
77: extern PetscErrorCode ComputeFunction(AppCtx*,Vec,Vec),FormInitialGuess(AppCtx*,Vec);
78: extern PetscErrorCode ComputeJacobian(AppCtx*,Vec,Mat);
80: int main(int argc,char **argv)
81: {
82: /* -------------- Data to define application problem ---------------- */
83: MPI_Comm comm; /* communicator */
84: KSP ksp; /* linear solver */
85: Vec X,Y,F; /* solution, update, residual vectors */
86: Mat J; /* Jacobian matrix */
87: AppCtx user; /* user-defined work context */
88: PetscInt Nx,Ny; /* number of preocessors in x- and y- directions */
89: PetscMPIInt size; /* number of processors */
90: PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
91: PetscInt m,N;
93: /* --------------- Data to define nonlinear solver -------------- */
94: PetscReal rtol = 1.e-8; /* relative convergence tolerance */
95: PetscReal xtol = 1.e-8; /* step convergence tolerance */
96: PetscReal ttol; /* convergence tolerance */
97: PetscReal fnorm,ynorm,xnorm; /* various vector norms */
98: PetscInt max_nonlin_its = 3; /* maximum number of iterations for nonlinear solver */
99: PetscInt max_functions = 50; /* maximum number of function evaluations */
100: PetscInt lin_its; /* number of linear solver iterations for each step */
101: PetscInt i; /* nonlinear solve iteration number */
102: PetscBool no_output = PETSC_FALSE; /* flag indicating whether to surpress output */
104: PetscInitialize(&argc,&argv,(char*)0,help);
105: comm = PETSC_COMM_WORLD;
106: PetscOptionsGetBool(NULL,NULL,"-no_output",&no_output,NULL);
108: /*
109: Initialize problem parameters
110: */
111: user.mx = 4; user.my = 4; user.param = 6.0;
113: PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,NULL);
114: PetscOptionsGetInt(NULL,NULL,"-my",&user.my,NULL);
115: PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL);
117: N = user.mx*user.my;
119: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: Create linear solver context
121: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123: KSPCreate(comm,&ksp);
125: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: Create vector data structures
127: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129: /*
130: Create distributed array (DMDA) to manage parallel grid and vectors
131: */
132: MPI_Comm_size(comm,&size);
133: Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
134: PetscOptionsGetInt(NULL,NULL,"-Nx",&Nx,NULL);
135: PetscOptionsGetInt(NULL,NULL,"-Ny",&Ny,NULL);
137: DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,&user.da);
138: DMSetFromOptions(user.da);
139: DMSetUp(user.da);
141: /*
142: Extract global and local vectors from DMDA; then duplicate for remaining
143: vectors that are the same types
144: */
145: DMCreateGlobalVector(user.da,&X);
146: DMCreateLocalVector(user.da,&user.localX);
147: VecDuplicate(X,&F);
148: VecDuplicate(X,&Y);
150: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151: Create matrix data structure for Jacobian
152: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153: /*
154: Note: For the parallel case, vectors and matrices MUST be partitioned
155: accordingly. When using distributed arrays (DMDAs) to create vectors,
156: the DMDAs determine the problem partitioning. We must explicitly
157: specify the local matrix dimensions upon its creation for compatibility
158: with the vector distribution. Thus, the generic MatCreate() routine
159: is NOT sufficient when working with distributed arrays.
161: Note: Here we only approximately preallocate storage space for the
162: Jacobian. See the users manual for a discussion of better techniques
163: for preallocating matrix memory.
164: */
165: if (size == 1) {
166: MatCreateSeqAIJ(comm,N,N,5,NULL,&J);
167: } else {
168: VecGetLocalSize(X,&m);
169: MatCreateAIJ(comm,m,m,N,N,5,NULL,3,NULL,&J);
170: }
172: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173: Customize linear solver; set runtime options
174: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
176: /*
177: Set runtime options (e.g.,-ksp_monitor -ksp_rtol <rtol> -ksp_type <type>)
178: */
179: KSPSetFromOptions(ksp);
181: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182: Evaluate initial guess
183: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185: FormInitialGuess(&user,X);
186: ComputeFunction(&user,X,F); /* Compute F(X) */
187: VecNorm(F,NORM_2,&fnorm); /* fnorm = || F || */
188: ttol = fnorm*rtol;
189: if (!no_output) PetscPrintf(comm,"Initial function norm = %g\n",(double)fnorm);
191: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192: Solve nonlinear system with a user-defined method
193: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195: /*
196: This solver is a very simplistic inexact Newton method, with no
197: no damping strategies or bells and whistles. The intent of this code
198: is merely to demonstrate the repeated solution with KSP of linear
199: systems with the same nonzero structure.
201: This is NOT the recommended approach for solving nonlinear problems
202: with PETSc! We urge users to employ the SNES component for solving
203: nonlinear problems whenever possible with application codes, as it
204: offers many advantages over coding nonlinear solvers independently.
205: */
207: for (i=0; i<max_nonlin_its; i++) {
209: /*
210: Compute the Jacobian matrix.
211: */
212: ComputeJacobian(&user,X,J);
214: /*
215: Solve J Y = F, where J is the Jacobian matrix.
216: - First, set the KSP linear operators. Here the matrix that
217: defines the linear system also serves as the preconditioning
218: matrix.
219: - Then solve the Newton system.
220: */
221: KSPSetOperators(ksp,J,J);
222: KSPSolve(ksp,F,Y);
223: KSPGetIterationNumber(ksp,&lin_its);
225: /*
226: Compute updated iterate
227: */
228: VecNorm(Y,NORM_2,&ynorm); /* ynorm = || Y || */
229: VecAYPX(Y,-1.0,X); /* Y <- X - Y */
230: VecCopy(Y,X); /* X <- Y */
231: VecNorm(X,NORM_2,&xnorm); /* xnorm = || X || */
232: if (!no_output) {
233: PetscPrintf(comm," linear solve iterations = %D, xnorm=%g, ynorm=%g\n",lin_its,(double)xnorm,(double)ynorm);
234: }
236: /*
237: Evaluate new nonlinear function
238: */
239: ComputeFunction(&user,X,F); /* Compute F(X) */
240: VecNorm(F,NORM_2,&fnorm); /* fnorm = || F || */
241: if (!no_output) {
242: PetscPrintf(comm,"Iteration %D, function norm = %g\n",i+1,(double)fnorm);
243: }
245: /*
246: Test for convergence
247: */
248: if (fnorm <= ttol) {
249: if (!no_output) {
250: PetscPrintf(comm,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)ttol);
251: }
252: break;
253: }
254: if (ynorm < xtol*(xnorm)) {
255: if (!no_output) {
256: PetscPrintf(comm,"Converged due to small update length: %g < %g * %g\n",(double)ynorm,(double)xtol,(double)xnorm);
257: }
258: break;
259: }
260: if (i > max_functions) {
261: if (!no_output) {
262: PetscPrintf(comm,"Exceeded maximum number of function evaluations: %D > %D\n",i,max_functions);
263: }
264: break;
265: }
266: }
267: PetscPrintf(comm,"Number of nonlinear iterations = %D\n",i);
269: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
270: Free work space. All PETSc objects should be destroyed when they
271: are no longer needed.
272: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
274: MatDestroy(&J)); PetscCall(VecDestroy(&Y);
275: VecDestroy(&user.localX)); PetscCall(VecDestroy(&X);
276: VecDestroy(&F);
277: KSPDestroy(&ksp)); PetscCall(DMDestroy(&user.da);
278: PetscFinalize();
279: return 0;
280: }
281: /* ------------------------------------------------------------------- */
282: /*
283: FormInitialGuess - Forms initial approximation.
285: Input Parameters:
286: user - user-defined application context
287: X - vector
289: Output Parameter:
290: X - vector
291: */
292: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
293: {
294: PetscInt i,j,row,mx,my,xs,ys,xm,ym,gxm,gym,gxs,gys;
295: PetscReal one = 1.0,lambda,temp1,temp,hx,hy;
296: PetscScalar *x;
298: mx = user->mx; my = user->my; lambda = user->param;
299: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
300: temp1 = lambda/(lambda + one);
302: /*
303: Get a pointer to vector data.
304: - For default PETSc vectors, VecGetArray() returns a pointer to
305: the data array. Otherwise, the routine is implementation dependent.
306: - You MUST call VecRestoreArray() when you no longer need access to
307: the array.
308: */
309: VecGetArray(X,&x);
311: /*
312: Get local grid boundaries (for 2-dimensional DMDA):
313: xs, ys - starting grid indices (no ghost points)
314: xm, ym - widths of local grid (no ghost points)
315: gxs, gys - starting grid indices (including ghost points)
316: gxm, gym - widths of local grid (including ghost points)
317: */
318: DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
319: DMDAGetGhostCorners(user->da,&gxs,&gys,NULL,&gxm,&gym,NULL);
321: /*
322: Compute initial guess over the locally owned part of the grid
323: */
324: for (j=ys; j<ys+ym; j++) {
325: temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
326: for (i=xs; i<xs+xm; i++) {
327: row = i - gxs + (j - gys)*gxm;
328: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
329: x[row] = 0.0;
330: continue;
331: }
332: x[row] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
333: }
334: }
336: /*
337: Restore vector
338: */
339: VecRestoreArray(X,&x);
340: return 0;
341: }
342: /* ------------------------------------------------------------------- */
343: /*
344: ComputeFunction - Evaluates nonlinear function, F(x).
346: Input Parameters:
347: . X - input vector
348: . user - user-defined application context
350: Output Parameter:
351: . F - function vector
352: */
353: PetscErrorCode ComputeFunction(AppCtx *user,Vec X,Vec F)
354: {
355: PetscInt i,j,row,mx,my,xs,ys,xm,ym,gxs,gys,gxm,gym;
356: PetscReal two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx,sc;
357: PetscScalar u,uxx,uyy,*x,*f;
358: Vec localX = user->localX;
360: mx = user->mx; my = user->my; lambda = user->param;
361: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
362: sc = hx*hy*lambda; hxdhy = hx/hy; hydhx = hy/hx;
364: /*
365: Scatter ghost points to local vector, using the 2-step process
366: DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
367: By placing code between these two statements, computations can be
368: done while messages are in transition.
369: */
370: DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
371: DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
373: /*
374: Get pointers to vector data
375: */
376: VecGetArray(localX,&x);
377: VecGetArray(F,&f);
379: /*
380: Get local grid boundaries
381: */
382: DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
383: DMDAGetGhostCorners(user->da,&gxs,&gys,NULL,&gxm,&gym,NULL);
385: /*
386: Compute function over the locally owned part of the grid
387: */
388: for (j=ys; j<ys+ym; j++) {
389: row = (j - gys)*gxm + xs - gxs - 1;
390: for (i=xs; i<xs+xm; i++) {
391: row++;
392: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
393: f[row] = x[row];
394: continue;
395: }
396: u = x[row];
397: uxx = (two*u - x[row-1] - x[row+1])*hydhx;
398: uyy = (two*u - x[row-gxm] - x[row+gxm])*hxdhy;
399: f[row] = uxx + uyy - sc*PetscExpScalar(u);
400: }
401: }
403: /*
404: Restore vectors
405: */
406: VecRestoreArray(localX,&x);
407: VecRestoreArray(F,&f);
408: PetscLogFlops(11.0*ym*xm);
409: return 0;
410: }
411: /* ------------------------------------------------------------------- */
412: /*
413: ComputeJacobian - Evaluates Jacobian matrix.
415: Input Parameters:
416: . x - input vector
417: . user - user-defined application context
419: Output Parameters:
420: . jac - Jacobian matrix
421: . flag - flag indicating matrix structure
423: Notes:
424: Due to grid point reordering with DMDAs, we must always work
425: with the local grid points, and then transform them to the new
426: global numbering with the "ltog" mapping
427: We cannot work directly with the global numbers for the original
428: uniprocessor grid!
429: */
430: PetscErrorCode ComputeJacobian(AppCtx *user,Vec X,Mat jac)
431: {
432: Vec localX = user->localX; /* local vector */
433: const PetscInt *ltog; /* local-to-global mapping */
434: PetscInt i,j,row,mx,my,col[5];
435: PetscInt xs,ys,xm,ym,gxs,gys,gxm,gym,grow;
436: PetscScalar two = 2.0,one = 1.0,lambda,v[5],hx,hy,hxdhy,hydhx,sc,*x;
437: ISLocalToGlobalMapping ltogm;
439: mx = user->mx; my = user->my; lambda = user->param;
440: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
441: sc = hx*hy; hxdhy = hx/hy; hydhx = hy/hx;
443: /*
444: Scatter ghost points to local vector, using the 2-step process
445: DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
446: By placing code between these two statements, computations can be
447: done while messages are in transition.
448: */
449: DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
450: DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
452: /*
453: Get pointer to vector data
454: */
455: VecGetArray(localX,&x);
457: /*
458: Get local grid boundaries
459: */
460: DMDAGetCorners(user->da,&xs,&ys,NULL,&xm,&ym,NULL);
461: DMDAGetGhostCorners(user->da,&gxs,&gys,NULL,&gxm,&gym,NULL);
463: /*
464: Get the global node numbers for all local nodes, including ghost points
465: */
466: DMGetLocalToGlobalMapping(user->da,<ogm);
467: ISLocalToGlobalMappingGetIndices(ltogm,<og);
469: /*
470: Compute entries for the locally owned part of the Jacobian.
471: - Currently, all PETSc parallel matrix formats are partitioned by
472: contiguous chunks of rows across the processors. The "grow"
473: parameter computed below specifies the global row number
474: corresponding to each local grid point.
475: - Each processor needs to insert only elements that it owns
476: locally (but any non-local elements will be sent to the
477: appropriate processor during matrix assembly).
478: - Always specify global row and columns of matrix entries.
479: - Here, we set all entries for a particular row at once.
480: */
481: for (j=ys; j<ys+ym; j++) {
482: row = (j - gys)*gxm + xs - gxs - 1;
483: for (i=xs; i<xs+xm; i++) {
484: row++;
485: grow = ltog[row];
486: /* boundary points */
487: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
488: MatSetValues(jac,1,&grow,1,&grow,&one,INSERT_VALUES);
489: continue;
490: }
491: /* interior grid points */
492: v[0] = -hxdhy; col[0] = ltog[row - gxm];
493: v[1] = -hydhx; col[1] = ltog[row - 1];
494: v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = grow;
495: v[3] = -hydhx; col[3] = ltog[row + 1];
496: v[4] = -hxdhy; col[4] = ltog[row + gxm];
497: MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
498: }
499: }
500: ISLocalToGlobalMappingRestoreIndices(ltogm,<og);
502: /*
503: Assemble matrix, using the 2-step process:
504: MatAssemblyBegin(), MatAssemblyEnd().
505: By placing code between these two statements, computations can be
506: done while messages are in transition.
507: */
508: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
509: VecRestoreArray(localX,&x);
510: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
512: return 0;
513: }
515: /*TEST
517: test:
519: TEST*/