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: const 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,PETSC_DECIDE,PETSC_DECIDE,N,N,&J);
102: MatSetFromOptions(J);
103: SNESSetJacobian(snes,J,J,FormJacobian1,&user);
104: }
106: /* Set PetscOptions, then solve nonlinear system */
107: SNESSetFromOptions(snes);
108: FormInitialGuess1(&user,x);
109: SNESSolve(snes,x);
110: SNESGetIterationNumber(snes,&its);
111: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %D\n",its);
113: /* Free data structures */
114: if (!matrix_free) {
115: MatDestroy(J);
116: }
117: VecDestroy(x);
118: VecDestroy(r);
119: VecDestroy(user.localX);
120: VecDestroy(user.localF);
121: SNESDestroy(snes);
122: DADestroy(user.da);
123: }
124: PetscFinalize();
126: return 0;
127: }/* -------------------- Form initial approximation ----------------- */
130: PetscErrorCode FormInitialGuess1(AppCtx *user,Vec X)
131: {
132: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xm,Ym,Xs,Ys;
134: PetscReal one = 1.0,lambda,temp1,temp,hx,hy;
135: PetscScalar *x;
136: Vec localX = user->localX;
138: mx = user->mx; my = user->my; lambda = user->param;
139: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
141: /* Get ghost points */
142: VecGetArray(localX,&x);
143: temp1 = lambda/(lambda + one);
144: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
145: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
147: /* Compute initial guess */
148: for (j=ys; j<ys+ym; j++) {
149: temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
150: for (i=xs; i<xs+xm; i++) {
151: row = i - Xs + (j - Ys)*Xm;
152: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
153: x[row] = 0.0;
154: continue;
155: }
156: x[row] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
157: }
158: }
159: VecRestoreArray(localX,&x);
161: /* Insert values into global vector */
162: DALocalToGlobal(user->da,localX,INSERT_VALUES,X);
163: return 0;
164: } /* -------------------- Evaluate Function F(x) --------------------- */
167: PetscErrorCode FormFunction1(SNES snes,Vec X,Vec F,void *ptr)
168: {
169: AppCtx *user = (AppCtx*)ptr;
171: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym;
172: PetscReal two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx,sc;
173: PetscScalar u,uxx,uyy,*x,*f;
174: Vec localX = user->localX,localF = user->localF;
176: mx = user->mx; my = user->my; lambda = user->param;
177: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
178: sc = hx*hy*lambda; hxdhy = hx/hy; hydhx = hy/hx;
180: /* Get ghost points */
181: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
182: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
183: VecGetArray(localX,&x);
184: VecGetArray(localF,&f);
185: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
186: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
188: /* Evaluate function */
189: for (j=ys; j<ys+ym; j++) {
190: row = (j - Ys)*Xm + xs - Xs - 1;
191: for (i=xs; i<xs+xm; i++) {
192: row++;
193: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
194: f[row] = x[row];
195: continue;
196: }
197: u = x[row];
198: uxx = (two*u - x[row-1] - x[row+1])*hydhx;
199: uyy = (two*u - x[row-Xm] - x[row+Xm])*hxdhy;
200: f[row] = uxx + uyy - sc*PetscExpScalar(u);
201: }
202: }
203: VecRestoreArray(localX,&x);
204: VecRestoreArray(localF,&f);
206: /* Insert values into global vector */
207: DALocalToGlobal(user->da,localF,INSERT_VALUES,F);
208: PetscLogFlops(11*ym*xm);
209: return 0;
210: } /* -------------------- Evaluate Jacobian F'(x) --------------------- */
213: PetscErrorCode FormJacobian1(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
214: {
215: AppCtx *user = (AppCtx*)ptr;
216: Mat jac = *J;
218: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym,col[5];
219: PetscInt nloc,*ltog,grow;
220: PetscScalar two = 2.0,one = 1.0,lambda,v[5],hx,hy,hxdhy,hydhx,sc,*x;
221: Vec localX = user->localX;
223: mx = user->mx; my = user->my; lambda = user->param;
224: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
225: sc = hx*hy; hxdhy = hx/hy; hydhx = hy/hx;
227: /* Get ghost points */
228: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
229: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
230: VecGetArray(localX,&x);
231: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
232: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
233: DAGetGlobalIndices(user->da,&nloc,<og);
235: /* Evaluate function */
236: for (j=ys; j<ys+ym; j++) {
237: row = (j - Ys)*Xm + xs - Xs - 1;
238: for (i=xs; i<xs+xm; i++) {
239: row++;
240: grow = ltog[row];
241: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
242: MatSetValues(jac,1,&grow,1,&grow,&one,INSERT_VALUES);
243: continue;
244: }
245: v[0] = -hxdhy; col[0] = ltog[row - Xm];
246: v[1] = -hydhx; col[1] = ltog[row - 1];
247: v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = grow;
248: v[3] = -hydhx; col[3] = ltog[row + 1];
249: v[4] = -hxdhy; col[4] = ltog[row + Xm];
250: MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
251: }
252: }
253: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
254: VecRestoreArray(X,&x);
255: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
256: *flag = SAME_NONZERO_PATTERN;
257: return 0;
258: }