Actual source code: ex1.c
1: /*$Id: ex1.c,v 1.82 2001/04/10 19:37:00 bsmith Exp $*/
3: /* Program usage: ex4 [-help] [all PETSc options] */
5: static char help[] = "Solves a nonlinear system on 1 processor with SNES. Wen
6: solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular domain.n
7: This example also illustrates the use of matrix coloring. Runtime options include:n
8: -par <parameter>, where <parameter> indicates the problem's nonlinearityn
9: problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)n
10: -mx <xg>, where <xg> = number of grid points in the x-directionn
11: -my <yg>, where <yg> = number of grid points in the y-directionnn";
13: /*T
14: Concepts: SNES^sequential Bratu example
15: Processors: 1
16: T*/
18: /* ------------------------------------------------------------------------
20: Solid Fuel Ignition (SFI) problem. This problem is modeled by
21: the partial differential equation
22:
23: -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
24:
25: with boundary conditions
26:
27: u = 0 for x = 0, x = 1, y = 0, y = 1.
28:
29: A finite difference approximation with the usual 5-point stencil
30: is used to discretize the boundary value problem to obtain a nonlinear
31: system of equations.
33: The parallel version of this code is snes/examples/tutorials/ex5.c
35: ------------------------------------------------------------------------- */
37: /*
38: Include "petscdraw.h" so that we can use PETSc drawing routines.
39: Include "petscsnes.h" so that we can use SNES solvers. Note that
40: this file automatically includes:
41: petsc.h - base PETSc routines petscvec.h - vectors
42: petscsys.h - system routines petscmat.h - matrices
43: petscis.h - index sets petscksp.h - Krylov subspace methods
44: petscviewer.h - viewers petscpc.h - preconditioners
45: petscsles.h - linear solvers
46: */
48: #include petscsnes.h
50: /*
51: User-defined application context - contains data needed by the
52: application-provided call-back routines, FormJacobian() and
53: FormFunction().
54: */
55: typedef struct {
56: double param; /* test problem parameter */
57: int mx; /* Discretization in x-direction */
58: int my; /* Discretization in y-direction */
59: } AppCtx;
61: /*
62: User-defined routines
63: */
64: extern int FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
65: extern int FormFunction(SNES,Vec,Vec,void*);
66: extern int FormInitialGuess(AppCtx*,Vec);
68: int main(int argc,char **argv)
69: {
70: SNES snes; /* nonlinear solver context */
71: Vec x,r; /* solution, residual vectors */
72: Mat J; /* Jacobian matrix */
73: AppCtx user; /* user-defined application context */
74: PetscDraw draw; /* drawing context */
75: int i,ierr,its,N,size,hist_its[50];
76: double bratu_lambda_max = 6.81,bratu_lambda_min = 0.,history[50];
77: MatFDColoring fdcoloring;
78: Scalar *array;
79: PetscTruth matrix_free,flg,fd_coloring;
81: PetscInitialize(&argc,&argv,(char *)0,help);
82: MPI_Comm_size(PETSC_COMM_WORLD,&size);
83: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
85: /*
86: Initialize problem parameters
87: */
88: user.mx = 4; user.my = 4; user.param = 6.0;
89: PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
90: PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
91: PetscOptionsGetDouble(PETSC_NULL,"-par",&user.param,PETSC_NULL);
92: if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
93: SETERRQ(1,"Lambda is out of range");
94: }
95: N = user.mx*user.my;
96:
97: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98: Create nonlinear solver context
99: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101: SNESCreate(PETSC_COMM_WORLD,SNES_NONLINEAR_EQUATIONS,&snes);
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Create vector data structures; set function evaluation routine
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,N,&x);
108: VecSetFromOptions(x);
109: VecDuplicate(x,&r);
111: /*
112: Set function evaluation routine and vector. Whenever the nonlinear
113: solver needs to evaluate the nonlinear function, it will call this
114: routine.
115: - Note that the final routine argument is the user-defined
116: context that provides application-specific data for the
117: function evaluation routine.
118: */
119: SNESSetFunction(snes,r,FormFunction,(void *)&user);
121: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122: Create matrix data structure; set Jacobian evaluation routine
123: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125: /*
126: Create matrix. Here we only approximately preallocate storage space
127: for the Jacobian. See the users manual for a discussion of better
128: techniques for preallocating matrix memory.
129: */
130: PetscOptionsHasName(PETSC_NULL,"-snes_mf",&matrix_free);
131: if (!matrix_free) {
132: MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,PETSC_NULL,&J);
133: }
135: /*
136: This option will cause the Jacobian to be computed via finite differences
137: efficiently using a coloring of the columns of the matrix.
138: */
139: PetscOptionsHasName(PETSC_NULL,"-snes_fd_coloring",&fd_coloring);
141: if (matrix_free && fd_coloring) SETERRQ(1,"Use only one of -snes_mf, -snes_fd_coloring options!n
142: You can do -snes_mf_operator -snes_fd_coloring");
144: if (fd_coloring) {
145: ISColoring iscoloring;
146: MatStructure str;
148: /*
149: This initializes the nonzero structure of the Jacobian. This is artificial
150: because clearly if we had a routine to compute the Jacobian we won't need
151: to use finite differences.
152: */
153: FormJacobian(snes,x,&J,&J,&str,&user);
155: /*
156: Color the matrix, i.e. determine groups of columns that share no common
157: rows. These columns in the Jacobian can all be computed simulataneously.
158: */
159: MatGetColoring(J,MATCOLORING_NATURAL,&iscoloring);
160: /*
161: Create the data structure that SNESDefaultComputeJacobianColor() uses
162: to compute the actual Jacobians via finite differences.
163: */
164: MatFDColoringCreate(J,iscoloring,&fdcoloring);
165: MatFDColoringSetFunction(fdcoloring,(int (*)(void))FormFunction,&user);
166: MatFDColoringSetFromOptions(fdcoloring);
167: /*
168: Tell SNES to use the routine SNESDefaultComputeJacobianColor()
169: to compute Jacobians.
170: */
171: SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,fdcoloring);
172: ISColoringDestroy(iscoloring);
173: }
174: /*
175: Set Jacobian matrix data structure and default Jacobian evaluation
176: routine. Whenever the nonlinear solver needs to compute the
177: Jacobian matrix, it will call this routine.
178: - Note that the final routine argument is the user-defined
179: context that provides application-specific data for the
180: Jacobian evaluation routine.
181: - The user can override with:
182: -snes_fd : default finite differencing approximation of Jacobian
183: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
184: (unless user explicitly sets preconditioner)
185: -snes_mf_operator : form preconditioning matrix as set by the user,
186: but use matrix-free approx for Jacobian-vector
187: products within Newton-Krylov method
188: */
189: else if (!matrix_free) {
190: SNESSetJacobian(snes,J,J,FormJacobian,(void*)&user);
191: }
193: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194: Customize nonlinear solver; set runtime options
195: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197: /*
198: Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
199: */
200: SNESSetFromOptions(snes);
202: /*
203: Set array that saves the function norms. This array is intended
204: when the user wants to save the convergence history for later use
205: rather than just to view the function norms via -snes_monitor.
206: */
207: SNESSetConvergenceHistory(snes,history,hist_its,50,PETSC_TRUE);
209: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210: Evaluate initial guess; then solve nonlinear system
211: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212: /*
213: Note: The user should initialize the vector, x, with the initial guess
214: for the nonlinear solver prior to calling SNESSolve(). In particular,
215: to employ an initial guess of zero, the user should explicitly set
216: this vector to zero by calling VecSet().
217: */
218: FormInitialGuess(&user,x);
219: SNESSolve(snes,x,&its);
220: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %dn",its);
222: /*
223: PetscDraw contour plot of solution
224: */
225: /* PetscDrawOpenX(PETSC_COMM_WORLD,0,"Solution",300,0,300,300,&draw); */
226: PetscDrawCreate(PETSC_COMM_WORLD,0,"Solution",300,0,300,300,&draw);
227: PetscDrawSetType(draw,PETSC_DRAW_X);
229: VecGetArray(x,&array);
230: PetscDrawTensorContour(draw,user.mx,user.my,0,0,array);
231: VecRestoreArray(x,&array);
233: /*
234: Print the convergence history. This is intended just to demonstrate
235: use of the data attained via SNESSetConvergenceHistory().
236: */
237: PetscOptionsHasName(PETSC_NULL,"-print_history",&flg);
238: if (flg) {
239: for (i=0; i<its+1; i++) {
240: PetscPrintf(PETSC_COMM_WORLD,"iteration %d: Linear iterations %d Function norm = %gn",i,hist_its[i],history[i]);
241: }
242: }
244: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
245: Free work space. All PETSc objects should be destroyed when they
246: are no longer needed.
247: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
249: if (!matrix_free) {
250: MatDestroy(J);
251: }
252: if (fd_coloring) {
253: MatFDColoringDestroy(fdcoloring);
254: }
255: VecDestroy(x);
256: VecDestroy(r);
257: PetscDrawDestroy(draw);
258: SNESDestroy(snes);
259: PetscFinalize();
261: return 0;
262: }
263: /* ------------------------------------------------------------------- */
264: /*
265: FormInitialGuess - Forms initial approximation.
267: Input Parameters:
268: user - user-defined application context
269: X - vector
271: Output Parameter:
272: X - vector
273: */
274: int FormInitialGuess(AppCtx *user,Vec X)
275: {
276: int i,j,row,mx,my,ierr;
277: double lambda,temp1,temp,hx,hy;
278: Scalar *x;
280: mx = user->mx;
281: my = user->my;
282: lambda = user->param;
284: hx = 1.0 / (double)(mx-1);
285: hy = 1.0 / (double)(my-1);
287: /*
288: Get a pointer to vector data.
289: - For default PETSc vectors, VecGetArray() returns a pointer to
290: the data array. Otherwise, the routine is implementation dependent.
291: - You MUST call VecRestoreArray() when you no longer need access to
292: the array.
293: */
294: VecGetArray(X,&x);
295: temp1 = lambda/(lambda + 1.0);
296: for (j=0; j<my; j++) {
297: temp = (double)(PetscMin(j,my-j-1))*hy;
298: for (i=0; i<mx; i++) {
299: row = i + j*mx;
300: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
301: x[row] = 0.0;
302: continue;
303: }
304: x[row] = temp1*sqrt(PetscMin((double)(PetscMin(i,mx-i-1))*hx,temp));
305: }
306: }
308: /*
309: Restore vector
310: */
311: VecRestoreArray(X,&x);
312: return 0;
313: }
314: /* ------------------------------------------------------------------- */
315: /*
316: FormFunction - Evaluates nonlinear function, F(x).
318: Input Parameters:
319: . snes - the SNES context
320: . X - input vector
321: . ptr - optional user-defined context, as set by SNESSetFunction()
323: Output Parameter:
324: . F - function vector
325: */
326: int FormFunction(SNES snes,Vec X,Vec F,void *ptr)
327: {
328: AppCtx *user = (AppCtx*)ptr;
329: int ierr,i,j,row,mx,my;
330: double two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx;
331: Scalar ut,ub,ul,ur,u,uxx,uyy,sc,*x,*f;
333: mx = user->mx;
334: my = user->my;
335: lambda = user->param;
336: hx = one / (double)(mx-1);
337: hy = one / (double)(my-1);
338: sc = hx*hy;
339: hxdhy = hx/hy;
340: hydhx = hy/hx;
342: /*
343: Get pointers to vector data
344: */
345: VecGetArray(X,&x);
346: VecGetArray(F,&f);
348: /*
349: Compute function
350: */
351: for (j=0; j<my; j++) {
352: for (i=0; i<mx; i++) {
353: row = i + j*mx;
354: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
355: f[row] = x[row];
356: continue;
357: }
358: u = x[row];
359: ub = x[row - mx];
360: ul = x[row - 1];
361: ut = x[row + mx];
362: ur = x[row + 1];
363: uxx = (-ur + two*u - ul)*hydhx;
364: uyy = (-ut + two*u - ub)*hxdhy;
365: f[row] = uxx + uyy - sc*lambda*PetscExpScalar(u);
366: }
367: }
369: /*
370: Restore vectors
371: */
372: VecRestoreArray(X,&x);
373: VecRestoreArray(F,&f);
374: return 0;
375: }
376: /* ------------------------------------------------------------------- */
377: /*
378: FormJacobian - Evaluates Jacobian matrix.
380: Input Parameters:
381: . snes - the SNES context
382: . x - input vector
383: . ptr - optional user-defined context, as set by SNESSetJacobian()
385: Output Parameters:
386: . A - Jacobian matrix
387: . B - optionally different preconditioning matrix
388: . flag - flag indicating matrix structure
389: */
390: int FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
391: {
392: AppCtx *user = (AppCtx*)ptr; /* user-defined applicatin context */
393: Mat jac = *J; /* Jacobian matrix */
394: int i,j,row,mx,my,col[5],ierr;
395: Scalar two = 2.0,one = 1.0,lambda,v[5],sc,*x;
396: double hx,hy,hxdhy,hydhx;
398: mx = user->mx;
399: my = user->my;
400: lambda = user->param;
401: hx = 1.0 / (double)(mx-1);
402: hy = 1.0 / (double)(my-1);
403: sc = hx*hy;
404: hxdhy = hx/hy;
405: hydhx = hy/hx;
407: /*
408: Get pointer to vector data
409: */
410: VecGetArray(X,&x);
412: /*
413: Compute entries of the Jacobian
414: */
415: for (j=0; j<my; j++) {
416: for (i=0; i<mx; i++) {
417: row = i + j*mx;
418: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
419: MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES);
420: continue;
421: }
422: v[0] = -hxdhy; col[0] = row - mx;
423: v[1] = -hydhx; col[1] = row - 1;
424: v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = row;
425: v[3] = -hydhx; col[3] = row + 1;
426: v[4] = -hxdhy; col[4] = row + mx;
427: MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES);
428: }
429: }
431: /*
432: Restore vector
433: */
434: VecRestoreArray(X,&x);
436: /*
437: Assemble matrix
438: */
439: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
440: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
442: /*
443: Set flag to indicate that the Jacobian matrix retains an identical
444: nonzero structure throughout all nonlinear iterations (although the
445: values of the entries change). Thus, we can save some work in setting
446: up the preconditioner (e.g., no need to redo symbolic factorization for
447: ILU/ICC preconditioners).
448: - If the nonzero structure of the matrix is different during
449: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
450: must be used instead. If you are unsure whether the matrix
451: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
452: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
453: believes your assertion and does not check the structure
454: of the matrix. If you erroneously claim that the structure
455: is the same when it actually is not, the new preconditioner
456: will not function correctly. Thus, use this optimization
457: feature with caution!
458: */
459: *flag = SAME_NONZERO_PATTERN;
460: return 0;
461: }