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.;
58: #if defined(PETSC_USE_LOG)
59: PetscLogStage stages[2];
60: #endif
62: PetscInitialize(&argc,&argv,(char *)0,help);
64: PetscLogStageRegister("stage 1",&stages[0]);
65: PetscLogStageRegister("stage 2",&stages[1]);
66: for (i=0; i<2; i++) {
67: PetscLogStagePush(stages[i]);
68: user.mx = 4; user.my = 4; user.param = 6.0;
69:
70: if (i!=0) {
71: PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
72: PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
73: PetscOptionsGetReal(PETSC_NULL,"-par",&user.param,PETSC_NULL);
74: if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
75: SETERRQ(1,"Lambda is out of range");
76: }
77: }
78: N = user.mx*user.my;
80: MPI_Comm_size(PETSC_COMM_WORLD,&size);
81: PetscOptionsGetInt(PETSC_NULL,"-Nx",&Nx,PETSC_NULL);
82: PetscOptionsGetInt(PETSC_NULL,"-Ny",&Ny,PETSC_NULL);
83: if (Nx*Ny != size && (Nx != PETSC_DECIDE || Ny != PETSC_DECIDE))
84: SETERRQ(1,"Incompatible number of processors: Nx * Ny != size");
85:
86: /* Set up distributed array */
87: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,user.mx,user.my,Nx,Ny,1,1,
88: PETSC_NULL,PETSC_NULL,&user.da);
89: DACreateGlobalVector(user.da,&x);
90: VecDuplicate(x,&r);
91: DACreateLocalVector(user.da,&user.localX);
92: VecDuplicate(user.localX,&user.localF);
94: /* Create nonlinear solver and set function evaluation routine */
95: SNESCreate(PETSC_COMM_WORLD,&snes);
96: SNESSetType(snes,type);
97: SNESSetFunction(snes,r,FormFunction1,&user);
99: /* Set default Jacobian evaluation routine. User can override with:
100: -snes_mf : matrix-free Newton-Krylov method with no preconditioning
101: (unless user explicitly sets preconditioner)
102: -snes_fd : default finite differencing approximation of Jacobian
103: */
104: PetscOptionsHasName(PETSC_NULL,"-snes_mf",&matrix_free);
105: if (!matrix_free) {
106: MatCreate(PETSC_COMM_WORLD,&J);
107: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,N,N);
108: MatSetFromOptions(J);
109: SNESSetJacobian(snes,J,J,FormJacobian1,&user);
110: }
112: /* Set PetscOptions, then solve nonlinear system */
113: SNESSetFromOptions(snes);
114: FormInitialGuess1(&user,x);
115: SNESSolve(snes,PETSC_NULL,x);
116: SNESGetIterationNumber(snes,&its);
117: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %D\n",its);
119: /* Free data structures */
120: if (!matrix_free) {
121: MatDestroy(J);
122: }
123: VecDestroy(x);
124: VecDestroy(r);
125: VecDestroy(user.localX);
126: VecDestroy(user.localF);
127: SNESDestroy(snes);
128: DADestroy(user.da);
129: }
130: PetscFinalize();
132: return 0;
133: }/* -------------------- Form initial approximation ----------------- */
136: PetscErrorCode FormInitialGuess1(AppCtx *user,Vec X)
137: {
138: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xm,Ym,Xs,Ys;
140: PetscReal one = 1.0,lambda,temp1,temp,hx,hy;
141: PetscScalar *x;
142: Vec localX = user->localX;
144: mx = user->mx; my = user->my; lambda = user->param;
145: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
147: /* Get ghost points */
148: VecGetArray(localX,&x);
149: temp1 = lambda/(lambda + one);
150: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
151: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
153: /* Compute initial guess */
154: for (j=ys; j<ys+ym; j++) {
155: temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
156: for (i=xs; i<xs+xm; i++) {
157: row = i - Xs + (j - Ys)*Xm;
158: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
159: x[row] = 0.0;
160: continue;
161: }
162: x[row] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
163: }
164: }
165: VecRestoreArray(localX,&x);
167: /* Insert values into global vector */
168: DALocalToGlobal(user->da,localX,INSERT_VALUES,X);
169: return 0;
170: } /* -------------------- Evaluate Function F(x) --------------------- */
173: PetscErrorCode FormFunction1(SNES snes,Vec X,Vec F,void *ptr)
174: {
175: AppCtx *user = (AppCtx*)ptr;
177: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym;
178: PetscReal two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx,sc;
179: PetscScalar u,uxx,uyy,*x,*f;
180: Vec localX = user->localX,localF = user->localF;
182: mx = user->mx; my = user->my; lambda = user->param;
183: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
184: sc = hx*hy*lambda; hxdhy = hx/hy; hydhx = hy/hx;
186: /* Get ghost points */
187: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
188: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
189: VecGetArray(localX,&x);
190: VecGetArray(localF,&f);
191: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
192: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
194: /* Evaluate function */
195: for (j=ys; j<ys+ym; j++) {
196: row = (j - Ys)*Xm + xs - Xs - 1;
197: for (i=xs; i<xs+xm; i++) {
198: row++;
199: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
200: f[row] = x[row];
201: continue;
202: }
203: u = x[row];
204: uxx = (two*u - x[row-1] - x[row+1])*hydhx;
205: uyy = (two*u - x[row-Xm] - x[row+Xm])*hxdhy;
206: f[row] = uxx + uyy - sc*PetscExpScalar(u);
207: }
208: }
209: VecRestoreArray(localX,&x);
210: VecRestoreArray(localF,&f);
212: /* Insert values into global vector */
213: DALocalToGlobal(user->da,localF,INSERT_VALUES,F);
214: PetscLogFlops(11*ym*xm);
215: return 0;
216: } /* -------------------- Evaluate Jacobian F'(x) --------------------- */
219: PetscErrorCode FormJacobian1(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
220: {
221: AppCtx *user = (AppCtx*)ptr;
222: Mat jac = *J;
224: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym,col[5];
225: PetscInt nloc,*ltog,grow;
226: PetscScalar two = 2.0,one = 1.0,lambda,v[5],hx,hy,hxdhy,hydhx,sc,*x;
227: Vec localX = user->localX;
229: mx = user->mx; my = user->my; lambda = user->param;
230: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
231: sc = hx*hy; hxdhy = hx/hy; hydhx = hy/hx;
233: /* Get ghost points */
234: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
235: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
236: VecGetArray(localX,&x);
237: DAGetCorners(user->da,&xs,&ys,0,&xm,&ym,0);
238: DAGetGhostCorners(user->da,&Xs,&Ys,0,&Xm,&Ym,0);
239: DAGetGlobalIndices(user->da,&nloc,<og);
241: /* Evaluate function */
242: for (j=ys; j<ys+ym; j++) {
243: row = (j - Ys)*Xm + xs - Xs - 1;
244: for (i=xs; i<xs+xm; i++) {
245: row++;
246: grow = ltog[row];
247: if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
248: MatSetValues(jac,1,&grow,1,&grow,&one,INSERT_VALUES);
249: continue;
250: }
251: v[0] = -hxdhy; col[0] = ltog[row - Xm];
252: v[1] = -hydhx; col[1] = ltog[row - 1];
253: v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = grow;
254: v[3] = -hydhx; col[3] = ltog[row + 1];
255: v[4] = -hxdhy; col[4] = ltog[row + Xm];
256: MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
257: }
258: }
259: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
260: VecRestoreArray(X,&x);
261: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
262: *flag = SAME_NONZERO_PATTERN;
263: return 0;
264: }