Actual source code: ex1.c
2: static char help[] = "Solves the nonlinear system, the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular domain.\n\
3: This example also illustrates the use of matrix coloring. Runtime options include:\n\
4: -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
5: problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\
6: -mx <xg>, where <xg> = number of grid points in the x-direction\n\
7: -my <yg>, where <yg> = number of grid points in the y-direction\n\n";
9: /*
11: Solid Fuel Ignition (SFI) problem. This problem is modeled by
12: the partial differential equation
14: -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
16: with boundary conditions
18: u = 0 for x = 0, x = 1, y = 0, y = 1.
20: A finite difference approximation with the usual 5-point stencil
21: is used to discretize the boundary value problem to obtain a nonlinear
22: system of equations.
24: The parallel version of this code is snes/tutorials/ex5.c
26: */
28: /*
29: Include "petscsnes.h" so that we can use SNES solvers. Note that
30: this file automatically includes:
31: petscsys.h - base PETSc routines petscvec.h - vectors
32: petscmat.h - matrices
33: petscis.h - index sets petscksp.h - Krylov subspace methods
34: petscviewer.h - viewers petscpc.h - preconditioners
35: petscksp.h - linear solvers
36: */
38: #include <petscsnes.h>
40: /*
41: User-defined application context - contains data needed by the
42: application-provided call-back routines, FormJacobian() and
43: FormFunction().
44: */
45: typedef struct {
46: PetscReal param; /* test problem parameter */
47: PetscInt mx; /* Discretization in x-direction */
48: PetscInt my; /* Discretization in y-direction */
49: } AppCtx;
51: /*
52: User-defined routines
53: */
54: extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
55: extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
56: extern PetscErrorCode FormInitialGuess(AppCtx *, Vec);
57: extern PetscErrorCode ConvergenceTest(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
58: extern PetscErrorCode ConvergenceDestroy(void *);
59: extern PetscErrorCode postcheck(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *);
61: int main(int argc, char **argv)
62: {
63: SNES snes; /* nonlinear solver context */
64: Vec x, r; /* solution, residual vectors */
65: Mat J; /* Jacobian matrix */
66: AppCtx user; /* user-defined application context */
67: PetscInt i, its, N, hist_its[50];
68: PetscMPIInt size;
69: PetscReal bratu_lambda_max = 6.81, bratu_lambda_min = 0., history[50];
70: MatFDColoring fdcoloring;
71: PetscBool matrix_free = PETSC_FALSE, flg, fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE, pc = PETSC_FALSE, prunejacobian = PETSC_FALSE;
72: KSP ksp;
73: PetscInt *testarray;
75: PetscFunctionBeginUser;
76: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
77: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
78: PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
80: /*
81: Initialize problem parameters
82: */
83: user.mx = 4;
84: user.my = 4;
85: user.param = 6.0;
86: PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, NULL));
87: PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, NULL));
88: PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL));
89: PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc", &pc, NULL));
90: PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda is out of range");
91: N = user.mx * user.my;
92: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_convergence_test", &use_convergence_test, NULL));
93: PetscCall(PetscOptionsGetBool(NULL, NULL, "-prune_jacobian", &prunejacobian, NULL));
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Create nonlinear solver context
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
101: if (pc) {
102: PetscCall(SNESSetType(snes, SNESNEWTONTR));
103: PetscCall(SNESNewtonTRSetPostCheck(snes, postcheck, NULL));
104: }
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Create vector data structures; set function evaluation routine
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110: PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
111: PetscCall(VecSetSizes(x, PETSC_DECIDE, N));
112: PetscCall(VecSetFromOptions(x));
113: PetscCall(VecDuplicate(x, &r));
115: /*
116: Set function evaluation routine and vector. Whenever the nonlinear
117: solver needs to evaluate the nonlinear function, it will call this
118: routine.
119: - Note that the final routine argument is the user-defined
120: context that provides application-specific data for the
121: function evaluation routine.
122: */
123: PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)&user));
125: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: Create matrix data structure; set Jacobian evaluation routine
127: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129: /*
130: Create matrix. Here we only approximately preallocate storage space
131: for the Jacobian. See the users manual for a discussion of better
132: techniques for preallocating matrix memory.
133: */
134: PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf", &matrix_free, NULL));
135: if (!matrix_free) {
136: PetscBool matrix_free_operator = PETSC_FALSE;
137: PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf_operator", &matrix_free_operator, NULL));
138: if (matrix_free_operator) matrix_free = PETSC_FALSE;
139: }
140: if (!matrix_free) PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, N, N, 5, NULL, &J));
142: /*
143: This option will cause the Jacobian to be computed via finite differences
144: efficiently using a coloring of the columns of the matrix.
145: */
146: PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_fd_coloring", &fd_coloring, NULL));
147: PetscCheck(!matrix_free || !fd_coloring, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Use only one of -snes_mf, -snes_fd_coloring options!\nYou can do -snes_mf_operator -snes_fd_coloring");
149: if (fd_coloring) {
150: ISColoring iscoloring;
151: MatColoring mc;
152: if (prunejacobian) {
153: /* Initialize x with random nonzero values so that the nonzeros in the Jacobian
154: can better reflect the sparsity structure of the Jacobian. */
155: PetscRandom rctx;
156: PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx));
157: PetscCall(PetscRandomSetInterval(rctx, 1.0, 2.0));
158: PetscCall(VecSetRandom(x, rctx));
159: PetscCall(PetscRandomDestroy(&rctx));
160: }
162: /*
163: This initializes the nonzero structure of the Jacobian. This is artificial
164: because clearly if we had a routine to compute the Jacobian we won't need
165: to use finite differences.
166: */
167: PetscCall(FormJacobian(snes, x, J, J, &user));
169: /*
170: Color the matrix, i.e. determine groups of columns that share no common
171: rows. These columns in the Jacobian can all be computed simultaneously.
172: */
173: PetscCall(MatColoringCreate(J, &mc));
174: PetscCall(MatColoringSetType(mc, MATCOLORINGSL));
175: PetscCall(MatColoringSetFromOptions(mc));
176: PetscCall(MatColoringApply(mc, &iscoloring));
177: PetscCall(MatColoringDestroy(&mc));
178: /*
179: Create the data structure that SNESComputeJacobianDefaultColor() uses
180: to compute the actual Jacobians via finite differences.
181: */
182: PetscCall(MatFDColoringCreate(J, iscoloring, &fdcoloring));
183: PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))FormFunction, &user));
184: PetscCall(MatFDColoringSetFromOptions(fdcoloring));
185: PetscCall(MatFDColoringSetUp(J, iscoloring, fdcoloring));
186: /*
187: Tell SNES to use the routine SNESComputeJacobianDefaultColor()
188: to compute Jacobians.
189: */
190: PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, fdcoloring));
191: PetscCall(ISColoringDestroy(&iscoloring));
192: if (prunejacobian) PetscCall(SNESPruneJacobianColor(snes, J, J));
193: }
194: /*
195: Set Jacobian matrix data structure and default Jacobian evaluation
196: routine. Whenever the nonlinear solver needs to compute the
197: Jacobian matrix, it will call this routine.
198: - Note that the final routine argument is the user-defined
199: context that provides application-specific data for the
200: Jacobian evaluation routine.
201: - The user can override with:
202: -snes_fd : default finite differencing approximation of Jacobian
203: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
204: (unless user explicitly sets preconditioner)
205: -snes_mf_operator : form preconditioning matrix as set by the user,
206: but use matrix-free approx for Jacobian-vector
207: products within Newton-Krylov method
208: */
209: else if (!matrix_free) {
210: PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, (void *)&user));
211: }
213: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214: Customize nonlinear solver; set runtime options
215: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217: /*
218: Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
219: */
220: PetscCall(SNESSetFromOptions(snes));
222: /*
223: Set array that saves the function norms. This array is intended
224: when the user wants to save the convergence history for later use
225: rather than just to view the function norms via -snes_monitor.
226: */
227: PetscCall(SNESSetConvergenceHistory(snes, history, hist_its, 50, PETSC_TRUE));
229: /*
230: Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the
231: user provided test before the specialized test. The convergence context is just an array to
232: test that it gets properly freed at the end
233: */
234: if (use_convergence_test) {
235: PetscCall(SNESGetKSP(snes, &ksp));
236: PetscCall(PetscMalloc1(5, &testarray));
237: PetscCall(KSPSetConvergenceTest(ksp, ConvergenceTest, testarray, ConvergenceDestroy));
238: }
240: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241: Evaluate initial guess; then solve nonlinear system
242: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243: /*
244: Note: The user should initialize the vector, x, with the initial guess
245: for the nonlinear solver prior to calling SNESSolve(). In particular,
246: to employ an initial guess of zero, the user should explicitly set
247: this vector to zero by calling VecSet().
248: */
249: PetscCall(FormInitialGuess(&user, x));
250: PetscCall(SNESSolve(snes, NULL, x));
251: PetscCall(SNESGetIterationNumber(snes, &its));
252: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
254: /*
255: Print the convergence history. This is intended just to demonstrate
256: use of the data attained via SNESSetConvergenceHistory().
257: */
258: PetscCall(PetscOptionsHasName(NULL, NULL, "-print_history", &flg));
259: if (flg) {
260: for (i = 0; i < its + 1; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iteration %" PetscInt_FMT ": Linear iterations %" PetscInt_FMT " Function norm = %g\n", i, hist_its[i], (double)history[i]));
261: }
263: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
264: Free work space. All PETSc objects should be destroyed when they
265: are no longer needed.
266: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
268: if (!matrix_free) PetscCall(MatDestroy(&J));
269: if (fd_coloring) PetscCall(MatFDColoringDestroy(&fdcoloring));
270: PetscCall(VecDestroy(&x));
271: PetscCall(VecDestroy(&r));
272: PetscCall(SNESDestroy(&snes));
273: PetscCall(PetscFinalize());
274: return 0;
275: }
277: /*
278: FormInitialGuess - Forms initial approximation.
280: Input Parameters:
281: user - user-defined application context
282: X - vector
284: Output Parameter:
285: X - vector
286: */
287: PetscErrorCode FormInitialGuess(AppCtx *user, Vec X)
288: {
289: PetscInt i, j, row, mx, my;
290: PetscReal lambda, temp1, temp, hx, hy;
291: PetscScalar *x;
293: PetscFunctionBeginUser;
294: mx = user->mx;
295: my = user->my;
296: lambda = user->param;
298: hx = 1.0 / (PetscReal)(mx - 1);
299: hy = 1.0 / (PetscReal)(my - 1);
301: /*
302: Get a pointer to vector data.
303: - For default PETSc vectors, VecGetArray() returns a pointer to
304: the data array. Otherwise, the routine is implementation dependent.
305: - You MUST call VecRestoreArray() when you no longer need access to
306: the array.
307: */
308: PetscCall(VecGetArray(X, &x));
309: temp1 = lambda / (lambda + 1.0);
310: for (j = 0; j < my; j++) {
311: temp = (PetscReal)(PetscMin(j, my - j - 1)) * hy;
312: for (i = 0; i < mx; i++) {
313: row = i + j * mx;
314: if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
315: x[row] = 0.0;
316: continue;
317: }
318: x[row] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, mx - i - 1)) * hx, temp));
319: }
320: }
322: /*
323: Restore vector
324: */
325: PetscCall(VecRestoreArray(X, &x));
326: PetscFunctionReturn(PETSC_SUCCESS);
327: }
329: /*
330: FormFunction - Evaluates nonlinear function, F(x).
332: Input Parameters:
333: . snes - the SNES context
334: . X - input vector
335: . ptr - optional user-defined context, as set by SNESSetFunction()
337: Output Parameter:
338: . F - function vector
339: */
340: PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr)
341: {
342: AppCtx *user = (AppCtx *)ptr;
343: PetscInt i, j, row, mx, my;
344: PetscReal two = 2.0, one = 1.0, lambda, hx, hy, hxdhy, hydhx;
345: PetscScalar ut, ub, ul, ur, u, uxx, uyy, sc, *f;
346: const PetscScalar *x;
348: PetscFunctionBeginUser;
349: mx = user->mx;
350: my = user->my;
351: lambda = user->param;
352: hx = one / (PetscReal)(mx - 1);
353: hy = one / (PetscReal)(my - 1);
354: sc = hx * hy;
355: hxdhy = hx / hy;
356: hydhx = hy / hx;
358: /*
359: Get pointers to vector data
360: */
361: PetscCall(VecGetArrayRead(X, &x));
362: PetscCall(VecGetArray(F, &f));
364: /*
365: Compute function
366: */
367: for (j = 0; j < my; j++) {
368: for (i = 0; i < mx; i++) {
369: row = i + j * mx;
370: if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
371: f[row] = x[row];
372: continue;
373: }
374: u = x[row];
375: ub = x[row - mx];
376: ul = x[row - 1];
377: ut = x[row + mx];
378: ur = x[row + 1];
379: uxx = (-ur + two * u - ul) * hydhx;
380: uyy = (-ut + two * u - ub) * hxdhy;
381: f[row] = uxx + uyy - sc * lambda * PetscExpScalar(u);
382: }
383: }
385: /*
386: Restore vectors
387: */
388: PetscCall(VecRestoreArrayRead(X, &x));
389: PetscCall(VecRestoreArray(F, &f));
390: PetscFunctionReturn(PETSC_SUCCESS);
391: }
393: /*
394: FormJacobian - Evaluates Jacobian matrix.
396: Input Parameters:
397: . snes - the SNES context
398: . x - input vector
399: . ptr - optional user-defined context, as set by SNESSetJacobian()
401: Output Parameters:
402: . A - Jacobian matrix
403: . B - optionally different preconditioning matrix
404: . flag - flag indicating matrix structure
405: */
406: PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr)
407: {
408: AppCtx *user = (AppCtx *)ptr; /* user-defined application context */
409: PetscInt i, j, row, mx, my, col[5];
410: PetscScalar two = 2.0, one = 1.0, lambda, v[5], sc;
411: const PetscScalar *x;
412: PetscReal hx, hy, hxdhy, hydhx;
414: PetscFunctionBeginUser;
415: mx = user->mx;
416: my = user->my;
417: lambda = user->param;
418: hx = 1.0 / (PetscReal)(mx - 1);
419: hy = 1.0 / (PetscReal)(my - 1);
420: sc = hx * hy;
421: hxdhy = hx / hy;
422: hydhx = hy / hx;
424: /*
425: Get pointer to vector data
426: */
427: PetscCall(VecGetArrayRead(X, &x));
429: /*
430: Compute entries of the Jacobian
431: */
432: for (j = 0; j < my; j++) {
433: for (i = 0; i < mx; i++) {
434: row = i + j * mx;
435: if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
436: PetscCall(MatSetValues(jac, 1, &row, 1, &row, &one, INSERT_VALUES));
437: continue;
438: }
439: v[0] = -hxdhy;
440: col[0] = row - mx;
441: v[1] = -hydhx;
442: col[1] = row - 1;
443: v[2] = two * (hydhx + hxdhy) - sc * lambda * PetscExpScalar(x[row]);
444: col[2] = row;
445: v[3] = -hydhx;
446: col[3] = row + 1;
447: v[4] = -hxdhy;
448: col[4] = row + mx;
449: PetscCall(MatSetValues(jac, 1, &row, 5, col, v, INSERT_VALUES));
450: }
451: }
453: /*
454: Restore vector
455: */
456: PetscCall(VecRestoreArrayRead(X, &x));
458: /*
459: Assemble matrix
460: */
461: PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
462: PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
464: if (jac != J) {
465: PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
466: PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
467: }
469: PetscFunctionReturn(PETSC_SUCCESS);
470: }
472: PetscErrorCode ConvergenceTest(KSP ksp, PetscInt it, PetscReal nrm, KSPConvergedReason *reason, void *ctx)
473: {
474: PetscFunctionBeginUser;
475: *reason = KSP_CONVERGED_ITERATING;
476: if (it > 1) {
477: *reason = KSP_CONVERGED_ITS;
478: PetscCall(PetscInfo(NULL, "User provided convergence test returning after 2 iterations\n"));
479: }
480: PetscFunctionReturn(PETSC_SUCCESS);
481: }
483: PetscErrorCode ConvergenceDestroy(void *ctx)
484: {
485: PetscFunctionBeginUser;
486: PetscCall(PetscInfo(NULL, "User provided convergence destroy called\n"));
487: PetscCall(PetscFree(ctx));
488: PetscFunctionReturn(PETSC_SUCCESS);
489: }
491: PetscErrorCode postcheck(SNES snes, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, void *ctx)
492: {
493: PetscReal norm;
494: Vec tmp;
496: PetscFunctionBeginUser;
497: PetscCall(VecDuplicate(x, &tmp));
498: PetscCall(VecWAXPY(tmp, -1.0, x, w));
499: PetscCall(VecNorm(tmp, NORM_2, &norm));
500: PetscCall(VecDestroy(&tmp));
501: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of search step %g\n", (double)norm));
502: PetscFunctionReturn(PETSC_SUCCESS);
503: }
505: /*TEST
507: build:
508: requires: !single
510: test:
511: args: -ksp_gmres_cgs_refinement_type refine_always
513: test:
514: suffix: 2
515: args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always
517: test:
518: suffix: 2a
519: filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault"
520: args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info
521: requires: defined(PETSC_USE_INFO)
523: test:
524: suffix: 2b
525: filter: grep -i "User provided convergence test" > /dev/null && echo "Found User provided convergence test"
526: args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info
527: requires: defined(PETSC_USE_INFO)
529: test:
530: suffix: 3
531: args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always
533: test:
534: suffix: 4
535: args: -pc -par 6.807 -snes_monitor -snes_converged_reason
537: test:
538: suffix: 5
539: args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always -prune_jacobian
540: output_file: output/ex1_3.out
542: test:
543: suffix: 6
544: args: -snes_monitor draw:image:testfile -viewer_view
546: TEST*/