Actual source code: ex22.c
2: static char help[] = "Solves PDE optimization problem.\n\n";
4: #include petscda.h
5: #include petscpf.h
6: #include petscsnes.h
7: #include petscdmmg.h
9: /*
11: w - design variables (what we change to get an optimal solution)
12: u - state variables (i.e. the PDE solution)
13: lambda - the Lagrange multipliers
15: U = (w u lambda)
17: fu, fw, flambda contain the gradient of L(w,u,lambda)
19: FU = (fw fu flambda)
21: In this example the PDE is
22: Uxx = 2,
23: u(0) = w(0), thus this is the free parameter
24: u(1) = 0
25: the function we wish to minimize is
26: \integral u^{2}
28: The exact solution for u is given by u(x) = x*x - 1.25*x + .25
30: Use the usual centered finite differences.
32: Note we treat the problem as non-linear though it happens to be linear
34: See ex21.c for the same code, but that does NOT interlaces the u and the lambda
36: The vectors u_lambda and fu_lambda contain the u and the lambda interlaced
37: */
39: typedef struct {
40: PetscViewer u_lambda_viewer;
41: PetscViewer fu_lambda_viewer;
42: } UserCtx;
50: int main(int argc,char **argv)
51: {
53: UserCtx user;
54: DA da;
55: DMMG *dmmg;
56: VecPack packer;
58: PetscInitialize(&argc,&argv,PETSC_NULL,help);
60: /* Hardwire several options; can be changed at command line */
61: PetscOptionsSetValue("-dmmg_grid_sequence",PETSC_NULL);
62: PetscOptionsSetValue("-ksp_type","fgmres");
63: PetscOptionsSetValue("-ksp_max_it","5");
64: PetscOptionsSetValue("-pc_mg_type","full");
65: PetscOptionsSetValue("-mg_coarse_ksp_type","gmres");
66: PetscOptionsSetValue("-mg_levels_ksp_type","gmres");
67: PetscOptionsSetValue("-mg_coarse_ksp_max_it","6");
68: PetscOptionsSetValue("-mg_levels_ksp_max_it","3");
69: PetscOptionsSetValue("-snes_mf_type","wp");
70: PetscOptionsSetValue("-snes_mf_compute_normu","no");
71: PetscOptionsSetValue("-snes_ls","basic");
72: PetscOptionsSetValue("-dmmg_jacobian_mf_fd",0);
73: /* PetscOptionsSetValue("-snes_ls","basicnonorms"); */
74: PetscOptionsInsert(&argc,&argv,PETSC_NULL);
76: /* Create a global vector that includes a single redundant array and two da arrays */
77: VecPackCreate(PETSC_COMM_WORLD,&packer);
78: VecPackAddArray(packer,1);
79: DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,-5,2,1,PETSC_NULL,&da);
80: VecPackAddDA(packer,da);
82: /* create graphics windows */
83: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u_lambda - state variables and Lagrange multipliers",-1,-1,-1,-1,&user.u_lambda_viewer);
84: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu_lambda - derivate w.r.t. state variables and Lagrange multipliers",-1,-1,-1,-1,&user.fu_lambda_viewer);
86: /* create nonlinear multi-level solver */
87: DMMGCreate(PETSC_COMM_WORLD,2,&user,&dmmg);
88: DMMGSetDM(dmmg,(DM)packer);
89: DMMGSetSNES(dmmg,FormFunction,PETSC_NULL);
90: /*
91: for (i=0; i<DMMGGetLevels(dmmg); i++) {
92: SNESSetMonitor(dmmg[i]->snes,Monitor,dmmg[i],0);
93: }*/
94: DMMGSolve(dmmg);
95: DMMGDestroy(dmmg);
97: DADestroy(da);
98: VecPackDestroy(packer);
99: PetscViewerDestroy(user.u_lambda_viewer);
100: PetscViewerDestroy(user.fu_lambda_viewer);
102: PetscFinalize();
103: return 0;
104: }
106: typedef struct {
107: PetscScalar u;
108: PetscScalar lambda;
109: } ULambda;
110:
111: /*
112: Evaluates FU = Gradiant(L(w,u,lambda))
114: This local function acts on the ghosted version of U (accessed via VecPackGetLocalVectors() and
115: VecPackScatter()) BUT the global, nonghosted version of FU (via VecPackGetAccess()).
117: */
118: PetscErrorCode FormFunction(SNES snes,Vec U,Vec FU,void* dummy)
119: {
120: DMMG dmmg = (DMMG)dummy;
122: PetscInt xs,xm,i,N,nredundant;
123: ULambda *u_lambda,*fu_lambda;
124: PetscScalar d,h,*w,*fw;
125: Vec vu_lambda,vfu_lambda;
126: DA da;
127: VecPack packer = (VecPack)dmmg->dm;
130: VecPackGetEntries(packer,&nredundant,&da);
131: VecPackGetLocalVectors(packer,&w,&vu_lambda);
132: VecPackScatter(packer,U,w,vu_lambda);
133: VecPackGetAccess(packer,FU,&fw,&vfu_lambda);
135: DAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
136: DAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0);
137: DAVecGetArray(da,vu_lambda,&u_lambda);
138: DAVecGetArray(da,vfu_lambda,&fu_lambda);
139: d = N-1.0;
140: h = 1.0/d;
142: /* derivative of L() w.r.t. w */
143: if (xs == 0) { /* only first processor computes this */
144: fw[0] = -2.0*d*u_lambda[0].lambda;
145: }
147: /* derivative of L() w.r.t. u */
148: for (i=xs; i<xs+xm; i++) {
149: if (i == 0) fu_lambda[0].lambda = h*u_lambda[0].u + 2.*d*u_lambda[0].lambda - d*u_lambda[1].lambda;
150: else if (i == 1) fu_lambda[1].lambda = 2.*h*u_lambda[1].u + 2.*d*u_lambda[1].lambda - d*u_lambda[2].lambda;
151: else if (i == N-1) fu_lambda[N-1].lambda = h*u_lambda[N-1].u + 2.*d*u_lambda[N-1].lambda - d*u_lambda[N-2].lambda;
152: else if (i == N-2) fu_lambda[N-2].lambda = 2.*h*u_lambda[N-2].u + 2.*d*u_lambda[N-2].lambda - d*u_lambda[N-3].lambda;
153: else fu_lambda[i].lambda = 2.*h*u_lambda[i].u - d*(u_lambda[i+1].lambda - 2.0*u_lambda[i].lambda + u_lambda[i-1].lambda);
154: }
156: /* derivative of L() w.r.t. lambda */
157: for (i=xs; i<xs+xm; i++) {
158: if (i == 0) fu_lambda[0].u = 2.0*d*(u_lambda[0].u - w[0]);
159: else if (i == N-1) fu_lambda[N-1].u = 2.0*d*u_lambda[N-1].u;
160: else fu_lambda[i].u = -(d*(u_lambda[i+1].u - 2.0*u_lambda[i].u + u_lambda[i-1].u) - 2.0*h);
161: }
163: DAVecRestoreArray(da,vu_lambda,&u_lambda);
164: DAVecRestoreArray(da,vfu_lambda,&fu_lambda);
165: VecPackRestoreLocalVectors(packer,&w,&vu_lambda);
166: VecPackRestoreAccess(packer,FU,&fw,&vfu_lambda);
167: PetscLogFlops(13*N);
168: return(0);
169: }
171: /*
172: Computes the exact solution
173: */
174: PetscErrorCode u_solution(void *dummy,PetscInt n,PetscScalar *x,PetscScalar *u)
175: {
176: PetscInt i;
178: for (i=0; i<n; i++) {
179: u[2*i] = x[i]*x[i] - 1.25*x[i] + .25;
180: }
181: return(0);
182: }
184: PetscErrorCode ExactSolution(VecPack packer,Vec U)
185: {
186: PF pf;
187: Vec x,u_global;
188: PetscScalar *w;
189: DA da;
191: PetscInt m;
194: VecPackGetEntries(packer,&m,&da);
196: PFCreate(PETSC_COMM_WORLD,1,1,&pf);
197: PFSetType(pf,PFQUICK,(void*)u_solution);
198: DAGetCoordinates(da,&x);
199: if (!x) {
200: DASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);
201: DAGetCoordinates(da,&x);
202: }
203: VecPackGetAccess(packer,U,&w,&u_global,0);
204: if (w) w[0] = .25;
205: PFApplyVec(pf,x,u_global);
206: PFDestroy(pf);
207: VecPackRestoreAccess(packer,U,&w,&u_global,0);
208: return(0);
209: }
212: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy)
213: {
214: DMMG dmmg = (DMMG)dummy;
215: UserCtx *user = (UserCtx*)dmmg->user;
217: PetscInt m,N;
218: PetscScalar *w,*dw;
219: Vec u_lambda,U,F,Uexact;
220: VecPack packer = (VecPack)dmmg->dm;
221: PetscReal norm;
222: DA da;
225: SNESGetSolution(snes,&U);
226: VecPackGetAccess(packer,U,&w,&u_lambda);
227: VecView(u_lambda,user->u_lambda_viewer);
228: VecPackRestoreAccess(packer,U,&w,&u_lambda);
230: SNESGetFunction(snes,&F,0,0);
231: VecPackGetAccess(packer,F,&w,&u_lambda);
232: /* VecView(u_lambda,user->fu_lambda_viewer); */
233: VecPackRestoreAccess(packer,U,&w,&u_lambda);
235: VecPackGetEntries(packer,&m,&da);
236: DAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0);
237: VecDuplicate(U,&Uexact);
238: ExactSolution(packer,Uexact);
239: VecAXPY(Uexact,-1.0,U);
240: VecPackGetAccess(packer,Uexact,&dw,&u_lambda);
241: VecStrideNorm(u_lambda,0,NORM_2,&norm);
242: norm = norm/sqrt(N-1.);
243: if (dw) PetscPrintf(dmmg->comm,"Norm of error %G Error at x = 0 %G\n",norm,PetscRealPart(dw[0]));
244: VecView(u_lambda,user->fu_lambda_viewer);
245: VecPackRestoreAccess(packer,Uexact,&dw,&u_lambda);
246: VecDestroy(Uexact);
247: return(0);
248: }