Actual source code: ex13.c
1: /*$Id: ex13.c,v 1.33 2001/08/07 21:31:12 bsmith Exp $*/
3: static char help[] = "This program is a replica of ex6.c except that it does 2 solves to avoid paging.n
4: This program demonstrates use of the SNES package to solve systems ofn
5: nonlinear equations in parallel, using 2-dimensional distributed arrays.n
6: The 2-dim Bratu (SFI - solid fuel ignition) test problem is used, wheren
7: analytic formation of the Jacobian is the default. The command linen
8: options are:n
9: -par <parameter>, where <parameter> indicates the problem's nonlinearityn
10: problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)n
11: -mx <xg>, where <xg> = number of grid points in the x-directionn
12: -my <yg>, where <yg> = number of grid points in the y-directionn
13: -Nx <npx>, where <npx> = number of processors in the x-directionn
14: -Ny <npy>, where <npy> = number of processors in the y-directionnn";
16: /*
17: 1) Solid Fuel Ignition (SFI) problem. This problem is modeled by
18: the partial differential equation
19:
20: -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
21:
22: with boundary conditions
23:
24: u = 0 for x = 0, x = 1, y = 0, y = 1.
25:
26: A finite difference approximation with the usual 5-point stencil
27: is used to discretize the boundary value problem to obtain a nonlinear
28: system of equations.
29: */
31: #include petscsnes.h
32: #include petscda.h
34: /* User-defined application context */
35: typedef struct {
36: PetscReal param; /* test problem parameter */
37: int mx,my; /* discretization in x, y directions */
38: Vec localX,localF; /* ghosted local vector */
39: DA da; /* distributed array data structure */
40: } AppCtx;
42: extern int FormFunction1(SNES,Vec,Vec,void*),FormInitialGuess1(AppCtx*,Vec);
43: extern int FormJacobian1(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
45: int main(int argc,char **argv)
46: {
47: SNES snes; /* nonlinear solver */
48: SNESType type = SNESLS; /* nonlinear solution method */
49: Vec x,r; /* solution, residual vectors */
50: Mat J; /* Jacobian matrix */
51: AppCtx user; /* user-defined work context */
52: int i,ierr,its,N,Nx = PETSC_DECIDE,Ny = PETSC_DECIDE;
53: PetscTruth matrix_free;
54: int size;
55: PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
57: PetscInitialize(&argc,&argv,(char *)0,help);
59: for (i=0; i<2; i++) {
60: PetscLogStagePush(i);
61: user.mx = 4; user.my = 4; user.param = 6.0;
62:
63: if (i!=0) {
64: PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
65: PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
66: PetscOptionsGetReal(PETSC_NULL,"-par",&user.param,PETSC_NULL);
67: if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
68: SETERRQ(1,"Lambda is out of range");
69: }
70: }
71: N = user.mx*user.my;
73: MPI_Comm_size(PETSC_COMM_WORLD,&size);
74: PetscOptionsGetInt(PETSC_NULL,"-Nx",&Nx,PETSC_NULL);
75: PetscOptionsGetInt(PETSC_NULL,"-Ny",&Ny,PETSC_NULL);
76: if (Nx*Ny != size && (Nx != PETSC_DECIDE || Ny != PETSC_DECIDE))
77: SETERRQ(1,"Incompatible number of processors: Nx * Ny != size");
78:
79: /* Set up distributed array */
80: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,user.mx,user.my,Nx,Ny,1,1,
81: PETSC_NULL,PETSC_NULL,&user.da);
82: DACreateGlobalVector(user.da,&x);
83: VecDuplicate(x,&r);
84: DACreateLocalVector(user.da,&user.localX);
85: VecDuplicate(user.localX,&user.localF);
87: /* Create nonlinear solver and set function evaluation routine */
88: SNESCreate(PETSC_COMM_WORLD,&snes);
89: SNESSetType(snes,type);
90: SNESSetFunction(snes,r,FormFunction1,&user);
92: /* Set default Jacobian evaluation routine. User can override with:
93: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
94: (unless user explicitly sets preconditioner)
95: -snes_fd : default finite differencing approximation of Jacobian
96: */
97: PetscOptionsHasName(PETSC_NULL,"-snes_mf",&matrix_free);
98: if (!matrix_free) {
99: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&J);
100: MatSetFromOptions(J);
101: SNESSetJacobian(snes,J,J,FormJacobian1,&user);
102: }
104: /* Set PetscOptions, then solve nonlinear system */
105: SNESSetFromOptions(snes);
106: FormInitialGuess1(&user,x);
107: SNESSolve(snes,x,&its);
108: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %dn",its);
110: /* Free data structures */
111: if (!matrix_free) {
112: MatDestroy(J);
113: }
114: VecDestroy(x);
115: VecDestroy(r);
116: VecDestroy(user.localX);
117: VecDestroy(user.localF);
118: SNESDestroy(snes);
119: DADestroy(user.da);
120: }
121: PetscFinalize();
123: return 0;
124: }/* -------------------- Form initial approximation ----------------- */
125: int FormInitialGuess1(AppCtx *user,Vec X)
126: {
127: int i,j,row,mx,my,ierr,xs,ys,xm,ym,Xm,Ym,Xs,Ys;
128: PetscReal one = 1.0,lambda,temp1,temp,hx,hy;
129: PetscScalar *x;
130: Vec localX = user->localX;
132: mx = user->mx; my = user->my; lambda = user->param;
133: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
135: /* Get ghost points */
136: VecGetArray(localX,&x);
137: temp1 = lambda/(lambda + one);
138: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
139: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
141: /* Compute initial guess */
142: for (j=ys; j<ys+ym; j++) {
143: temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
144: for (i=xs; i<xs+xm; i++) {
145: row = i - Xs + (j - Ys)*Xm;
146: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
147: x[row] = 0.0;
148: continue;
149: }
150: x[row] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
151: }
152: }
153: VecRestoreArray(localX,&x);
155: /* Insert values into global vector */
156: DALocalToGlobal(user->da,localX,INSERT_VALUES,X);
157: return 0;
158: } /* -------------------- Evaluate Function F(x) --------------------- */
159: int FormFunction1(SNES snes,Vec X,Vec F,void *ptr)
160: {
161: AppCtx *user = (AppCtx*)ptr;
162: int ierr,i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym;
163: PetscReal two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx,sc;
164: PetscScalar u,uxx,uyy,*x,*f;
165: Vec localX = user->localX,localF = user->localF;
167: mx = user->mx; my = user->my; lambda = user->param;
168: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
169: sc = hx*hy*lambda; hxdhy = hx/hy; hydhx = hy/hx;
171: /* Get ghost points */
172: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
173: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
174: VecGetArray(localX,&x);
175: VecGetArray(localF,&f);
176: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
177: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
179: /* Evaluate function */
180: for (j=ys; j<ys+ym; j++) {
181: row = (j - Ys)*Xm + xs - Xs - 1;
182: for (i=xs; i<xs+xm; i++) {
183: row++;
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: uxx = (two*u - x[row-1] - x[row+1])*hydhx;
190: uyy = (two*u - x[row-Xm] - x[row+Xm])*hxdhy;
191: f[row] = uxx + uyy - sc*PetscExpScalar(u);
192: }
193: }
194: VecRestoreArray(localX,&x);
195: VecRestoreArray(localF,&f);
197: /* Insert values into global vector */
198: DALocalToGlobal(user->da,localF,INSERT_VALUES,F);
199: PetscLogFlops(11*ym*xm);
200: return 0;
201: } /* -------------------- Evaluate Jacobian F'(x) --------------------- */
202: int FormJacobian1(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
203: {
204: AppCtx *user = (AppCtx*)ptr;
205: Mat jac = *J;
206: int ierr,i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym,col[5];
207: int nloc,*ltog,grow;
208: PetscScalar two = 2.0,one = 1.0,lambda,v[5],hx,hy,hxdhy,hydhx,sc,*x;
209: Vec localX = user->localX;
211: mx = user->mx; my = user->my; lambda = user->param;
212: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
213: sc = hx*hy; hxdhy = hx/hy; hydhx = hy/hx;
215: /* Get ghost points */
216: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
217: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
218: VecGetArray(localX,&x);
219: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
220: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
221: DAGetGlobalIndices(user->da,&nloc,<og);
223: /* Evaluate function */
224: for (j=ys; j<ys+ym; j++) {
225: row = (j - Ys)*Xm + xs - Xs - 1;
226: for (i=xs; i<xs+xm; i++) {
227: row++;
228: grow = ltog[row];
229: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
230: MatSetValues(jac,1,&grow,1,&grow,&one,INSERT_VALUES);
231: continue;
232: }
233: v[0] = -hxdhy; col[0] = ltog[row - Xm];
234: v[1] = -hydhx; col[1] = ltog[row - 1];
235: v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = grow;
236: v[3] = -hydhx; col[3] = ltog[row + 1];
237: v[4] = -hxdhy; col[4] = ltog[row + Xm];
238: MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
239: }
240: }
241: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
242: VecRestoreArray(X,&x);
243: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
244: *flag = SAME_NONZERO_PATTERN;
245: return 0;
246: }