Actual source code: ex4.c
1: /*$Id: ex4.c,v 1.64 2001/08/07 21:31:12 bsmith Exp $*/
3: /* NOTE: THIS PROGRAM HAS NOT YET BEEN SET UP IN TUTORIAL STYLE. */
5: static char help[] ="This program demonstrates use of the SNES package to solve systems of nonlinear equations on a single processor.n
6: Both of these examples employn
7: sparse storage of the Jacobian matrices and are taken from the MINPACK-2n
8: test suite. By default the Bratu (SFI - solid fuel ignition) test problemn
9: is solved. The command line options are:n
10: -cavity : Solve FDC (flow in a driven cavity) problemn
11: -par <parameter>, where <parameter> indicates the problem's nonlinearityn
12: problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)n
13: problem FDC: <parameter> = Reynolds number (par > 0)n
14: -mx <xg>, where <xg> = number of grid points in the x-directionn
15: -my <yg>, where <yg> = number of grid points in the y-directionnn";
17: /*
18: 1) 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.
30:
31: 2) Flow in a Driven Cavity (FDC) problem. The problem is
32: formulated as a boundary value problem, which is discretized by
33: standard finite difference approximations to obtain a system of
34: nonlinear equations.
35: */
37: #include petscsnes.h
39: typedef struct {
40: PetscReal param; /* test problem parameter */
41: int mx; /* Discretization in x-direction */
42: int my; /* Discretization in y-direction */
43: } AppCtx;
45: extern int FormJacobian1(SNES,Vec,Mat*,Mat*,MatStructure*,void*),
46: FormFunction1(SNES,Vec,Vec,void*),
47: FormInitialGuess1(AppCtx*,Vec);
48: extern int FormJacobian2(SNES,Vec,Mat*,Mat*,MatStructure*,void*),
49: FormFunction2(SNES,Vec,Vec,void*),
50: FormInitialGuess2(AppCtx*,Vec);
51:
52: int main(int argc, char **argv)
53: {
54: SNES snes; /* SNES context */
55: SNESType method = SNESLS; /* default nonlinear solution method */
56: Vec x, r; /* solution, residual vectors */
57: Mat J; /* Jacobian matrix */
58: AppCtx user; /* user-defined application context */
59: PetscDraw draw; /* drawing context */
60: int ierr, its, N, nfails;
61: PetscTruth flg,cavity;
62: PetscReal bratu_lambda_max = 6.81, bratu_lambda_min = 0.;
63: PetscScalar *xvalues;
65: PetscInitialize(&argc, &argv,(char *)0,help);
66: /* PetscDrawOpenX(PETSC_COMM_WORLD,0,"Solution",300,0,300,300,&draw); */
67: PetscDrawCreate(PETSC_COMM_WORLD,0,"Solution",300,0,300,300,&draw);
68: PetscDrawSetType(draw,PETSC_DRAW_X);
70: user.mx = 4;
71: user.my = 4;
72: PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
73: PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
74: PetscOptionsHasName(PETSC_NULL,"-cavity",&cavity);
75: if (cavity) user.param = 100.0;
76: else user.param = 6.0;
77: PetscOptionsGetReal(PETSC_NULL,"-par",&user.param,PETSC_NULL);
78: if (!cavity && (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min)) {
79: SETERRQ(1,"Lambda is out of range");
80: }
81: N = user.mx*user.my;
82:
83: /* Set up data structures */
84: VecCreateSeq(PETSC_COMM_SELF,N,&x);
85: VecDuplicate(x,&r);
86: MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,5,PETSC_NULL,&J);
88: /* Create nonlinear solver */
89: SNESCreate(PETSC_COMM_WORLD,&snes);
90: SNESSetType(snes,method);
92: /* Set various routines and compute initial guess for nonlinear solver */
93: if (cavity){
94: FormInitialGuess2(&user,x);
95: SNESSetFunction(snes,r,FormFunction2,(void *)&user);
96: SNESSetJacobian(snes,J,J,FormJacobian2,(void *)&user);
97: } else {
98: FormInitialGuess1(&user,x);
99: SNESSetFunction(snes,r,FormFunction1,(void *)&user);
100: SNESSetJacobian(snes,J,J,FormJacobian1,(void *)&user);
101: }
103: /* Set options and solve nonlinear system */
104: SNESSetFromOptions(snes);
105: SNESSolve(snes,x,&its);
106: SNESGetNumberUnsuccessfulSteps(snes,&nfails);
108: PetscPrintf(PETSC_COMM_SELF,"number of Newton iterations = %d, ",its);
109: PetscPrintf(PETSC_COMM_SELF,"number of unsuccessful steps = %dnn",nfails);
110: VecGetArray(x,&xvalues);
111: PetscDrawTensorContour(draw,user.mx,user.my,0,0,xvalues);
112: VecRestoreArray(x,&xvalues);
114: /* Free data structures */
115: VecDestroy(x); VecDestroy(r);
116: MatDestroy(J); SNESDestroy(snes);
117: PetscDrawDestroy(draw);
118: PetscFinalize();
120: return 0;
121: }
122: /* ------------------------------------------------------------------ */
123: /* Bratu (Solid Fuel Ignition) Test Problem */
124: /* ------------------------------------------------------------------ */
126: /* -------------------- Form initial approximation ----------------- */
128: int FormInitialGuess1(AppCtx *user,Vec X)
129: {
130: int i, j, row, mx, my, ierr;
131: PetscReal lambda, temp1, temp, hx, hy, hxdhy, hydhx,sc;
132: PetscScalar *x;
134: mx = user->mx;
135: my = user->my;
136: lambda = user->param;
138: hx = 1.0 / (PetscReal)(mx-1);
139: hy = 1.0 / (PetscReal)(my-1);
140: sc = hx*hy;
141: hxdhy = hx/hy;
142: hydhx = hy/hx;
144: VecGetArray(X,&x);
145: temp1 = lambda/(lambda + 1.0);
146: for (j=0; j<my; j++) {
147: temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
148: for (i=0; i<mx; i++) {
149: row = i + j*mx;
150: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
151: x[row] = 0.0;
152: continue;
153: }
154: x[row] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
155: }
156: }
157: VecRestoreArray(X,&x);
158: return 0;
159: }
160: /* -------------------- Evaluate Function F(x) --------------------- */
161:
162: int FormFunction1(SNES snes,Vec X,Vec F,void *ptr)
163: {
164: AppCtx *user = (AppCtx*)ptr;
165: int ierr, i, j, row, mx, my;
166: PetscReal two = 2.0, one = 1.0, lambda,hx, hy, hxdhy, hydhx;
167: PetscScalar ut, ub, ul, ur, u, uxx, uyy, sc,*x,*f;
169: mx = user->mx;
170: my = user->my;
171: lambda = user->param;
173: hx = one / (PetscReal)(mx-1);
174: hy = one / (PetscReal)(my-1);
175: sc = hx*hy;
176: hxdhy = hx/hy;
177: hydhx = hy/hx;
179: VecGetArray(X,&x);
180: VecGetArray(F,&f);
181: for (j=0; j<my; j++) {
182: for (i=0; i<mx; i++) {
183: row = i + j*mx;
184: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
185: f[row] = x[row];
186: continue;
187: }
188: u = x[row];
189: ub = x[row - mx];
190: ul = x[row - 1];
191: ut = x[row + mx];
192: ur = x[row + 1];
193: uxx = (-ur + two*u - ul)*hydhx;
194: uyy = (-ut + two*u - ub)*hxdhy;
195: f[row] = uxx + uyy - sc*lambda*exp(u);
196: }
197: }
198: VecRestoreArray(X,&x);
199: VecRestoreArray(F,&f);
200: return 0;
201: }
202: /* -------------------- Evaluate Jacobian F'(x) -------------------- */
204: int FormJacobian1(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
205: {
206: AppCtx *user = (AppCtx*)ptr;
207: Mat jac = *J;
208: int i, j, row, mx, my, col[5], ierr;
209: PetscScalar two = 2.0, one = 1.0, lambda, v[5],sc, *x;
210: PetscReal hx, hy, hxdhy, hydhx;
213: mx = user->mx;
214: my = user->my;
215: lambda = user->param;
217: hx = 1.0 / (PetscReal)(mx-1);
218: hy = 1.0 / (PetscReal)(my-1);
219: sc = hx*hy;
220: hxdhy = hx/hy;
221: hydhx = hy/hx;
223: VecGetArray(X,&x);
224: for (j=0; j<my; j++) {
225: for (i=0; i<mx; i++) {
226: row = i + j*mx;
227: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
228: MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES);
229: continue;
230: }
231: v[0] = -hxdhy; col[0] = row - mx;
232: v[1] = -hydhx; col[1] = row - 1;
233: v[2] = two*(hydhx + hxdhy) - sc*lambda*exp(x[row]); col[2] = row;
234: v[3] = -hydhx; col[3] = row + 1;
235: v[4] = -hxdhy; col[4] = row + mx;
236: MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES);
237: }
238: }
239: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
240: VecRestoreArray(X,&x);
241: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
242: *flag = SAME_NONZERO_PATTERN;
243: return 0;
244: }
245: /* ------------------------------------------------------------------ */
246: /* Driven Cavity Test Problem */
247: /* ------------------------------------------------------------------ */
249: /* -------------------- Form initial approximation ----------------- */
251: int FormInitialGuess2(AppCtx *user,Vec X)
252: {
253: int ierr, i, j, row, mx, my;
254: PetscScalar xx,yy,*x;
255: PetscReal hx, hy;
257: mx = user->mx;
258: my = user->my;
260: /* Test for invalid input parameters */
261: if ((mx <= 0) || (my <= 0)) SETERRQ(1,0);
263: hx = 1.0 / (PetscReal)(mx-1);
264: hy = 1.0 / (PetscReal)(my-1);
266: VecGetArray(X,&x);
267: yy = 0.0;
268: for (j=0; j<my; j++) {
269: xx = 0.0;
270: for (i=0; i<mx; i++) {
271: row = i + j*mx;
272: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
273: x[row] = 0.0;
274: }
275: else {
276: x[row] = - xx*(1.0 - xx)*yy*(1.0 - yy);
277: }
278: xx = xx + hx;
279: }
280: yy = yy + hy;
281: }
282: VecRestoreArray(X,&x);
283: return 0;
284: }
285: /* -------------------- Evaluate Function F(x) --------------------- */
287: int FormFunction2(SNES snes,Vec X,Vec F,void *pptr)
288: {
289: AppCtx *user = (AppCtx*)pptr;
290: int i, j, row, mx, my, ierr;
291: PetscScalar two = 2.0, zero = 0.0, pb, pbb,pbr, pl,pll,p,pr,prr;
292: PetscScalar ptl,pt,ptt,dpdy,dpdx,pblap,ptlap,rey,pbl,ptr,pllap,plap,prlap;
293: PetscScalar *x,*f, hx2, hy2, hxhy2;
294: PetscReal hx, hy;
296: mx = user->mx;
297: my = user->my;
298: hx = 1.0 / (PetscReal)(mx-1);
299: hy = 1.0 / (PetscReal)(my-1);
300: hx2 = hx*hx;
301: hy2 = hy*hy;
302: hxhy2 = hx2*hy2;
303: rey = user->param;
305: VecGetArray(X,&x);
306: VecGetArray(F,&f);
307: for (j=0; j<my; j++) {
308: for (i=0; i<mx; i++) {
309: row = i + j*mx;
310: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
311: f[row] = x[row];
312: continue;
313: }
314: if (i == 1 || j == 1) {
315: pbl = zero;
316: }
317: else {
318: pbl = x[row-mx-1];
319: }
320: if (j == 1) {
321: pb = zero;
322: pbb = x[row];
323: }
324: else if (j == 2) {
325: pb = x[row-mx];
326: pbb = zero;
327: }
328: else {
329: pb = x[row-mx];
330: pbb = x[row-2*mx];
331: }
332: if (j == 1 || i == mx-2) {
333: pbr = zero;
334: }
335: else {
336: pbr = x[row-mx+1];
337: }
338: if (i == 1) {
339: pl = zero;
340: pll = x[row];
341: }
342: else if (i == 2) {
343: pl = x[row-1];
344: pll = zero;
345: }
346: else {
347: pl = x[row-1];
348: pll = x[row-2];
349: }
350: p = x[row];
351: if (i == mx-3) {
352: pr = x[row+1];
353: prr = zero;
354: }
355: else if (i == mx-2) {
356: pr = zero;
357: prr = x[row];
358: }
359: else {
360: pr = x[row+1];
361: prr = x[row+2];
362: }
363: if (j == my-2 || i == 1) {
364: ptl = zero;
365: }
366: else {
367: ptl = x[row+mx-1];
368: }
369: if (j == my-3) {
370: pt = x[row+mx];
371: ptt = zero;
372: }
373: else if (j == my-2) {
374: pt = zero;
375: ptt = x[row] + two*hy;
376: }
377: else {
378: pt = x[row+mx];
379: ptt = x[row+2*mx];
380: }
381: if (j == my-2 || i == mx-2) {
382: ptr = zero;
383: }
384: else {
385: ptr = x[row+mx+1];
386: }
388: dpdy = (pt - pb)/(two*hy);
389: dpdx = (pr - pl)/(two*hx);
391: /* Laplacians at each point in the 5 point stencil */
392: pblap = (pbr - two*pb + pbl)/hx2 + (p - two*pb + pbb)/hy2;
393: pllap = (p - two*pl + pll)/hx2 + (ptl - two*pl + pbl)/hy2;
394: plap = (pr - two*p + pl)/hx2 + (pt - two*p + pb)/hy2;
395: prlap = (prr - two*pr + p)/hx2 + (ptr - two*pr + pbr)/hy2;
396: ptlap = (ptr - two*pt + ptl)/hx2 + (ptt - two*pt + p)/hy2;
398: f[row] = hxhy2*((prlap - two*plap + pllap)/hx2
399: + (ptlap - two*plap + pblap)/hy2
400: - rey*(dpdy*(prlap - pllap)/(two*hx)
401: - dpdx*(ptlap - pblap)/(two*hy)));
402: }
403: }
404: VecRestoreArray(X,&x);
405: VecRestoreArray(F,&f);
406: return 0;
407: }
408: /* -------------------- Evaluate Jacobian F'(x) -------------------- */
410: int FormJacobian2(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *pptr)
411: {
412: AppCtx *user = (AppCtx*)pptr;
413: int i, j, row, mx, my, col, ierr;
414: PetscScalar two = 2.0, one = 1.0, zero = 0.0, pb, pbb,pbr, pl,pll,p,pr,prr;
415: PetscScalar ptl,pt,ptt,dpdy,dpdx,pblap,ptlap,rey,pbl,ptr,pllap,plap,prlap;
416: PetscScalar val,four = 4.0, three = 3.0,*x;
417: PetscReal hx, hy,hx2, hy2, hxhy2;
419: mx = user->mx;
420: my = user->my;
421: hx = 1.0 / (PetscReal)(mx-1);
422: hy = 1.0 / (PetscReal)(my-1);
423: hx2 = hx*hx;
424: hy2 = hy*hy;
425: hxhy2 = hx2*hy2;
426: rey = user->param;
428: MatZeroEntries(*J);
429: VecGetArray(X,&x);
430: for (j=0; j<my; j++) {
431: for (i=0; i<mx; i++) {
432: row = i + j*mx;
433: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
434: MatSetValues(*J,1,&row,1,&row,&one,ADD_VALUES);
435: continue;
436: }
437: if (i == 1 || j == 1) {
438: pbl = zero;
439: }
440: else {
441: pbl = x[row-mx-1];
442: }
443: if (j == 1) {
444: pb = zero;
445: pbb = x[row];
446: }
447: else if (j == 2) {
448: pb = x[row-mx];
449: pbb = zero;
450: }
451: else {
452: pb = x[row-mx];
453: pbb = x[row-2*mx];
454: }
455: if (j == 1 || i == mx-2) {
456: pbr = zero;
457: }
458: else {
459: pbr = x[row-mx+1];
460: }
461: if (i == 1) {
462: pl = zero;
463: pll = x[row];
464: }
465: else if (i == 2) {
466: pl = x[row-1];
467: pll = zero;
468: }
469: else {
470: pl = x[row-1];
471: pll = x[row-2];
472: }
473: p = x[row];
474: if (i == mx-3) {
475: pr = x[row+1];
476: prr = zero;
477: }
478: else if (i == mx-2) {
479: pr = zero;
480: prr = x[row];
481: }
482: else {
483: pr = x[row+1];
484: prr = x[row+2];
485: }
486: if (j == my-2 || i == 1) {
487: ptl = zero;
488: }
489: else {
490: ptl = x[row+mx-1];
491: }
492: if (j == my-3) {
493: pt = x[row+mx];
494: ptt = zero;
495: }
496: else if (j == my-2) {
497: pt = zero;
498: ptt = x[row] + two*hy;
499: }
500: else {
501: pt = x[row+mx];
502: ptt = x[row+2*mx];
503: }
504: if (j == my-2 || i == mx-2) {
505: ptr = zero;
506: }
507: else {
508: ptr = x[row+mx+1];
509: }
511: dpdy = (pt - pb)/(two*hy);
512: dpdx = (pr - pl)/(two*hx);
514: /* Laplacians at each point in the 5 point stencil */
515: pblap = (pbr - two*pb + pbl)/hx2 + (p - two*pb + pbb)/hy2;
516: pllap = (p - two*pl + pll)/hx2 + (ptl - two*pl + pbl)/hy2;
517: plap = (pr - two*p + pl)/hx2 + (pt - two*p + pb)/hy2;
518: prlap = (prr - two*pr + p)/hx2 + (ptr - two*pr + pbr)/hy2;
519: ptlap = (ptr - two*pt + ptl)/hx2 + (ptt - two*pt + p)/hy2;
521: if (j > 2) {
522: val = hxhy2*(one/hy2/hy2 - rey*dpdx/hy2/(two*hy));
523: col = row - 2*mx;
524: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
525: }
526: if (i > 2) {
527: val = hxhy2*(one/hx2/hx2 + rey*dpdy/hx2/(two*hx));
528: col = row - 2;
529: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
530: }
531: if (i < mx-3) {
532: val = hxhy2*(one/hx2/hx2 - rey*dpdy/hx2/(two*hx));
533: col = row + 2;
534: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
535: }
536: if (j < my-3) {
537: val = hxhy2*(one/hy2/hy2 + rey*dpdx/hy2/(two*hy));
538: col = row + 2*mx;
539: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
540: }
541: if (i != 1 && j != 1) {
542: val = hxhy2*(two/hy2/hx2 + rey*(dpdy/hy2/(two*hx) - dpdx/hx2/(two*hy)));
543: col = row - mx - 1;
544: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
545: }
546: if (j != 1 && i != mx-2) {
547: val = hxhy2*(two/hy2/hx2 - rey*(dpdy/hy2/(two*hx) + dpdx/hx2/(two*hy)));
548: col = row - mx + 1;
549: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
550: }
551: if (j != my-2 && i != 1) {
552: val = hxhy2*(two/hy2/hx2 + rey*(dpdy/hy2/(two*hx) + dpdx/hx2/(two*hy)));
553: col = row + mx - 1;
554: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
555: }
556: if (j != my-2 && i != mx-2) {
557: val = hxhy2*(two/hy2/hx2 - rey*(dpdy/hy2/(two*hx) - dpdx/hx2/(two*hy)));
558: col = row + mx + 1;
559: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
560: }
561: if (j != 1) {
562: val = hxhy2*(-four*(one/hy2/hx2 + one/hy2/hy2)
563: + rey*((prlap - pllap)/(two*hx)/(two*hy)
564: + dpdx*(one/hx2 + one/hy2)/hy));
565: col = row - mx;
566: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
567: }
568: if (i != 1) {
569: val = hxhy2*(-four*(one/hy2/hx2 + one/hx2/hx2)
570: - rey*((ptlap - pblap)/(two*hx)/(two*hy)
571: + dpdy*(one/hx2 + one/hy2)/hx));
572: col = row - 1;
573: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
574: }
575: if (i != mx-2) {
576: val = hxhy2*(-four*(one/hy2/hx2 + one/hx2/hx2)
577: + rey*((ptlap - pblap)/(two*hx)/(two*hy)
578: + dpdy*(one/hx2 + one/hy2)/hx));
579: col = row + 1;
580: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
581: }
582: if (j != my-2) {
583: val = hxhy2*(-four*(one/hy2/hx2 + one/hy2/hy2)
584: - rey*((prlap - pllap)/(two*hx)/(two*hy)
585: + dpdx*(one/hx2 + one/hy2)/hy));
586: col = row + mx;
587: MatSetValues(*J,1,&row,1,&col,&val,ADD_VALUES);
588: }
589: val = hxhy2*(two*(four/hx2/hy2 + three/hx2/hx2 + three/hy2/hy2));
590: MatSetValues(*J,1,&row,1,&row,&val,ADD_VALUES);
591: if (j == 1) {
592: val = hxhy2*(one/hy2/hy2 - rey*(dpdx/hy2/(two*hy)));
593: MatSetValues(*J,1,&row,1,&row,&val,ADD_VALUES);
594: }
595: if (i == 1) {
596: val = hxhy2*(one/hx2/hx2 + rey*(dpdy/hx2/(two*hx)));
597: MatSetValues(*J,1,&row,1,&row,&val,ADD_VALUES);
598: }
599: if (i == mx-2) {
600: val = hxhy2*(one/hx2/hx2 - rey*(dpdy/hx2/(two*hx)));
601: MatSetValues(*J,1,&row,1,&row,&val,ADD_VALUES);
602: }
603: if (j == my-2) {
604: val = hxhy2*(one/hy2/hy2 + rey*(dpdx/hy2/(two*hy)));
605: MatSetValues(*J,1,&row,1,&row,&val,ADD_VALUES);
606: }
607: }
608: }
609: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
610: VecRestoreArray(X,&x);
611: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
612: *flag = SAME_NONZERO_PATTERN;
613: return 0;
614: }