Actual source code: ex25.c
1: /* $Id: ex18.c,v 1.23 2001/08/07 21:31:17 bsmith Exp $ */
4: static char help[] ="Minimum surface problemn
5: Uses 2-dimensional distributed arrays.n
6: n
7: Solves the linear systems via multilevel methods n
8: nn";
10: /*T
11: Concepts: SNES^solving a system of nonlinear equations
12: Concepts: DA^using distributed arrays
13: Concepts: multigrid;
14: Processors: n
15: T*/
17: /*
18:
19: This example models the partial differential equation
20:
21: - Div((1 + ||GRAD T||^2)^(1/2) (GRAD T)) = 0.
22:
23:
24: in the unit square, which is uniformly discretized in each of x and
25: y in this simple encoding. The degrees of freedom are vertex centered
26:
27: A finite difference approximation with the usual 5-point stencil
28: is used to discretize the boundary value problem to obtain a
29: nonlinear system of equations.
30:
31: */
33: #include petscsnes.h
34: #include petscda.h
35: #include petscmg.h
37: extern int FormFunction(SNES,Vec,Vec,void*);
38: extern int FormFunctionLocal(DALocalInfo*,PetscScalar**,PetscScalar**,void*);
40: int main(int argc,char **argv)
41: {
42: DMMG *dmmg;
43: SNES snes;
44: int ierr,its,lits;
45: PetscReal litspit;
46: DA da;
48: PetscInitialize(&argc,&argv,PETSC_NULL,help);
51: /*
52: Create the multilevel DA data structure
53: */
54: DMMGCreate(PETSC_COMM_WORLD,3,0,&dmmg);
56: /*
57: Set the DA (grid structure) for the grids.
58: */
59: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,-5,-5,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
60: DMMGSetDM(dmmg,(DM)da);
61: DADestroy(da);
63: /*
64: Process adiC: FormFunctionLocal FormFunctionLocali
66: Create the nonlinear solver, and tell the DMMG structure to use it
67: */
68: /* DMMGSetSNES(dmmg,FormFunction,0); */
69: DMMGSetSNESLocal(dmmg,FormFunctionLocal,0,ad_FormFunctionLocal,0);
71: /*
72: PreLoadBegin() means that the following section of code is run twice. The first time
73: through the flag PreLoading is on this the nonlinear solver is only run for a single step.
74: The second time through (the actually timed code) the maximum iterations is set to 10
75: Preload of the executable is done to eliminate from the timing the time spent bring the
76: executable into memory from disk (paging in).
77: */
78: PreLoadBegin(PETSC_TRUE,"Solve");
79: DMMGSolve(dmmg);
80: PreLoadEnd();
81: snes = DMMGGetSNES(dmmg);
82: SNESGetIterationNumber(snes,&its);
83: SNESGetNumberLinearIterations(snes,&lits);
84: litspit = ((PetscReal)lits)/((PetscReal)its);
85: PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %dn",its);
86: PetscPrintf(PETSC_COMM_WORLD,"Number of Linear iterations = %dn",lits);
87: PetscPrintf(PETSC_COMM_WORLD,"Average Linear its / Newton = %en",litspit);
89: DMMGDestroy(dmmg);
90: PetscFinalize();
92: return 0;
93: }
94: /* -------------------- Evaluate Function F(x) --------------------- */
95: int FormFunction(SNES snes,Vec T,Vec F,void* ptr)
96: {
97: DMMG dmmg = (DMMG)ptr;
98: int ierr,i,j,mx,my,xs,ys,xm,ym;
99: PetscScalar hx,hy;
100: PetscScalar **t,**f,gradup,graddown,gradleft,gradright,gradx,grady;
101: PetscScalar coeffup,coeffdown,coeffleft,coeffright;
102: Vec localT;
105: DAGetLocalVector((DA)dmmg->dm,&localT);
106: DAGetInfo((DA)dmmg->dm,PETSC_NULL,&mx,&my,0,0,0,0,0,0,0,0);
107: hx = 1.0/(PetscReal)(mx-1); hy = 1.0/(PetscReal)(my-1);
108:
109: /* Get ghost points */
110: DAGlobalToLocalBegin((DA)dmmg->dm,T,INSERT_VALUES,localT);
111: DAGlobalToLocalEnd((DA)dmmg->dm,T,INSERT_VALUES,localT);
112: DAGetCorners((DA)dmmg->dm,&xs,&ys,0,&xm,&ym,0);
113: DAVecGetArray((DA)dmmg->dm,localT,(void**)&t);
114: DAVecGetArray((DA)dmmg->dm,F,(void**)&f);
116: /* Evaluate function */
117: for (j=ys; j<ys+ym; j++) {
118: for (i=xs; i<xs+xm; i++) {
120: if (i == 0 || i == mx-1 || j == 0 || j == my-1) {
122: f[j][i] = t[j][i] - (1.0 - (2.0*hx*(PetscReal)i - 1.0)*(2.0*hx*(PetscReal)i - 1.0));
123:
124: } else {
126: gradup = (t[j+1][i] - t[j][i])/hy;
127: graddown = (t[j][i] - t[j-1][i])/hy;
128: gradright = (t[j][i+1] - t[j][i])/hx;
129: gradleft = (t[j][i] - t[j][i-1])/hx;
131: gradx = .5*(t[j][i+1] - t[j][i-1])/hx;
132: grady = .5*(t[j+1][i] - t[j-1][i])/hy;
134: coeffup = 1.0/PetscSqrtScalar(1.0 + gradup*gradup + gradx*gradx);
135: coeffdown = 1.0/PetscSqrtScalar(1.0 + graddown*graddown + gradx*gradx);
137: coeffleft = 1.0/PetscSqrtScalar(1.0 + gradleft*gradleft + grady*grady);
138: coeffright = 1.0/PetscSqrtScalar(1.0 + gradright*gradright + grady*grady);
140: f[j][i] = (coeffup*gradup - coeffdown*graddown)*hx + (coeffright*gradright - coeffleft*gradleft)*hy;
141:
142: }
144: }
145: }
146: DAVecRestoreArray((DA)dmmg->dm,localT,(void**)&t);
147: DAVecRestoreArray((DA)dmmg->dm,F,(void**)&f);
148: DARestoreLocalVector((DA)dmmg->dm,&localT);
149: return(0);
150: }
152: int FormFunctionLocal(DALocalInfo *info,PetscScalar **t,PetscScalar **f,void *ptr)
153: {
154: int i,j;
155: PetscScalar hx,hy;
156: PetscScalar gradup,graddown,gradleft,gradright,gradx,grady;
157: PetscScalar coeffup,coeffdown,coeffleft,coeffright;
160: hx = 1.0/(PetscReal)(info->mx-1); hy = 1.0/(PetscReal)(info->my-1);
161:
162: /* Evaluate function */
163: for (j=info->ys; j<info->ys+info->ym; j++) {
164: for (i=info->xs; i<info->xs+info->xm; i++) {
166: if (i == 0 || i == info->mx-1 || j == 0 || j == info->my-1) {
168: f[j][i] = t[j][i] - (1.0 - (2.0*hx*(PetscReal)i - 1.0)*(2.0*hx*(PetscReal)i - 1.0));
169:
170: } else {
172: gradup = (t[j+1][i] - t[j][i])/hy;
173: graddown = (t[j][i] - t[j-1][i])/hy;
174: gradright = (t[j][i+1] - t[j][i])/hx;
175: gradleft = (t[j][i] - t[j][i-1])/hx;
177: gradx = .5*(t[j][i+1] - t[j][i-1])/hx;
178: grady = .5*(t[j+1][i] - t[j-1][i])/hy;
180: coeffup = 1.0/PetscSqrtScalar(1.0 + gradup*gradup + gradx*gradx);
181: coeffdown = 1.0/PetscSqrtScalar(1.0 + graddown*graddown + gradx*gradx);
183: coeffleft = 1.0/PetscSqrtScalar(1.0 + gradleft*gradleft + grady*grady);
184: coeffright = 1.0/PetscSqrtScalar(1.0 + gradright*gradright + grady*grady);
186: f[j][i] = (coeffup*gradup - coeffdown*graddown)*hx + (coeffright*gradright - coeffleft*gradleft)*hy;
187:
188: }
190: }
191: }
192: return(0);
193: }