Actual source code: ex9.c
1: /*$Id: ex9.c,v 1.52 2001/08/07 21:31:12 bsmith Exp $*/
3: static char help[] = "This program demonstrates use of the SNES package. Solve systems ofn
4: nonlinear equations in parallel. This example uses matrix-free Krylovn
5: Newton methods with no preconditioner to solve the Bratu (SFI - solid fueln
6: ignition) test problem. The command line options are:n
7: -par <parameter>, where <parameter> indicates the problem's nonlinearityn
8: problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)n
9: -mx <xg>, where <xg> = number of grid points in the x-directionn
10: -my <yg>, where <yg> = number of grid points in the y-directionn
11: -mz <zg>, where <zg> = number of grid points in the z-directionnn";
13: /*
14: 1) Solid Fuel Ignition (SFI) problem. This problem is modeled by
15: the partial differential equation
16:
17: -Laplacian u - lambda*exp(u) = 0, 0 < x,y,z < 1,
18:
19: with boundary conditions
20:
21: u = 0 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1.
22:
23: A finite difference approximation with the usual 7-point stencil
24: is used to discretize the boundary value problem to obtain a nonlinear
25: system of equations.
26: */
28: #include petscsnes.h
29: #include petscda.h
31: typedef struct {
32: PetscReal param; /* test problem nonlinearity parameter */
33: int mx,my,mz; /* discretization in x,y,z-directions */
34: Vec localX,localF; /* ghosted local vectors */
35: DA da; /* distributed array datastructure */
36: } AppCtx;
38: extern int FormFunction1(SNES,Vec,Vec,void*),FormInitialGuess1(AppCtx*,Vec);
40: int main(int argc,char **argv)
41: {
42: SNES snes; /* nonlinear solver */
43: SLES sles; /* linear solver */
44: PC pc; /* preconditioner */
45: Mat J; /* Jacobian matrix */
46: AppCtx user; /* user-defined application context */
47: Vec x,r; /* vectors */
48: DAStencilType stencil = DA_STENCIL_BOX;
49: int ierr,its;
50: PetscTruth flg;
51: int Nx = PETSC_DECIDE,Ny = PETSC_DECIDE,Nz = PETSC_DECIDE;
52: PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
54: PetscInitialize(&argc,&argv,(char *)0,help);
55: PetscOptionsHasName(PETSC_NULL,"-star",&flg);
56: if (flg) stencil = DA_STENCIL_STAR;
58: user.mx = 4;
59: user.my = 4;
60: user.mz = 4;
61: user.param = 6.0;
62: PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
63: PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
64: PetscOptionsGetInt(PETSC_NULL,"-mz",&user.mz,PETSC_NULL);
65: PetscOptionsGetInt(PETSC_NULL,"-Nx",&Nx,PETSC_NULL);
66: PetscOptionsGetInt(PETSC_NULL,"-Ny",&Ny,PETSC_NULL);
67: PetscOptionsGetInt(PETSC_NULL,"-Nz",&Nz,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: /* Set up distributed array */
74: DACreate3d(PETSC_COMM_WORLD,DA_NONPERIODIC,stencil,user.mx,user.my,user.mz,
75: Nx,Ny,Nz,1,1,PETSC_NULL,PETSC_NULL,PETSC_NULL,&user.da);
76: DACreateGlobalVector(user.da,&x);
77: VecDuplicate(x,&r);
78: DACreateLocalVector(user.da,&user.localX);
79: VecDuplicate(user.localX,&user.localF);
81: /* Create nonlinear solver */
82: SNESCreate(PETSC_COMM_WORLD,&snes);
83: /* Set various routines and options */
84: SNESSetFunction(snes,r,FormFunction1,(void *)&user);
85: MatCreateSNESMF(snes,x,&J);
86: SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,&user);
87: SNESSetFromOptions(snes);
89: /* Force no preconditioning to be used. Note that this overrides whatever
90: choices may have been specified in the options database. */
91: SNESGetSLES(snes,&sles);
92: SLESGetPC(sles,&pc);
93: PCSetType(pc,PCNONE);
95: /* Solve nonlinear system */
96: FormInitialGuess1(&user,x);
97: SNESSolve(snes,x,&its);
98: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %dn",its);
100: /* Free data structures */
101: VecDestroy(user.localX);
102: VecDestroy(user.localF);
103: DADestroy(user.da);
104: VecDestroy(x); VecDestroy(r);
105: MatDestroy(J); SNESDestroy(snes);
107: PetscFinalize();
108: return 0;
109: }/* -------------------- Form initial approximation ----------------- */
110: int FormInitialGuess1(AppCtx *user,Vec X)
111: {
112: int i,j,k,loc,mx,my,mz,ierr,xs,ys,zs,xm,ym,zm,Xm,Ym,Zm,Xs,Ys,Zs,base1;
113: PetscReal one = 1.0,lambda,temp1,temp,Hx,Hy;
114: PetscScalar *x;
115: Vec localX = user->localX;
117: mx = user->mx; my = user->my; mz = user->mz; lambda = user->param;
118: Hx = one / (PetscReal)(mx-1); Hy = one / (PetscReal)(my-1);
120: ierr = VecGetArray(localX,&x);
121: temp1 = lambda/(lambda + one);
122: ierr = DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);
123: ierr = DAGetGhostCorners(user->da,&Xs,&Ys,&Zs,&Xm,&Ym,&Zm);
124:
125: for (k=zs; k<zs+zm; k++) {
126: base1 = (Xm*Ym)*(k-Zs);
127: for (j=ys; j<ys+ym; j++) {
128: temp = (PetscReal)(PetscMin(j,my-j-1))*Hy;
129: for (i=xs; i<xs+xm; i++) {
130: loc = base1 + i-Xs + (j-Ys)*Xm;
131: if (i == 0 || j == 0 || k == 0 || i==mx-1 || j==my-1 || k==mz-1) {
132: x[loc] = 0.0;
133: continue;
134: }
135: x[loc] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*Hx,temp));
136: }
137: }
138: }
140: VecRestoreArray(localX,&x);
141: /* stick values into global vector */
142: DALocalToGlobal(user->da,localX,INSERT_VALUES,X);
143: return 0;
144: }/* -------------------- Evaluate Function F(x) --------------------- */
145: int FormFunction1(SNES snes,Vec X,Vec F,void *ptr)
146: {
147: AppCtx *user = (AppCtx*)ptr;
148: int ierr,i,j,k,loc,mx,my,mz,xs,ys,zs,xm,ym,zm,Xs,Ys,Zs,Xm,Ym,Zm;
149: int base1,base2;
150: PetscReal two = 2.0,one = 1.0,lambda,Hx,Hy,Hz,HxHzdHy,HyHzdHx,HxHydHz;
151: PetscScalar u,uxx,uyy,sc,*x,*f,uzz;
152: Vec localX = user->localX,localF = user->localF;
154: mx = user->mx; my = user->my; mz = user->mz; lambda = user->param;
155: Hx = one / (PetscReal)(mx-1);
156: Hy = one / (PetscReal)(my-1);
157: Hz = one / (PetscReal)(mz-1);
158: sc = Hx*Hy*Hz*lambda; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
159: HxHydHz = Hx*Hy/Hz;
161: DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
162: DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
163: VecGetArray(localX,&x);
164: VecGetArray(localF,&f);
166: DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);
167: DAGetGhostCorners(user->da,&Xs,&Ys,&Zs,&Xm,&Ym,&Zm);
169: for (k=zs; k<zs+zm; k++) {
170: base1 = (Xm*Ym)*(k-Zs);
171: for (j=ys; j<ys+ym; j++) {
172: base2 = base1 + (j-Ys)*Xm;
173: for (i=xs; i<xs+xm; i++) {
174: loc = base2 + (i-Xs);
175: if (i == 0 || j == 0 || k== 0 || i == mx-1 || j == my-1 || k == mz-1) {
176: f[loc] = x[loc];
177: }
178: else {
179: u = x[loc];
180: uxx = (two*u - x[loc-1] - x[loc+1])*HyHzdHx;
181: uyy = (two*u - x[loc-Xm] - x[loc+Xm])*HxHzdHy;
182: uzz = (two*u - x[loc-Xm*Ym] - x[loc+Xm*Ym])*HxHydHz;
183: f[loc] = uxx + uyy + uzz - sc*PetscExpScalar(u);
184: }
185: }
186: }
187: }
188: VecRestoreArray(localX,&x);
189: VecRestoreArray(localF,&f);
190: /* stick values into global vector */
191: DALocalToGlobal(user->da,localF,INSERT_VALUES,F);
192: PetscLogFlops(11*ym*xm*zm);
193: return 0;
194: }
195:
200: