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