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