Actual source code: ex24.c
2: static char help[] = "Solves PDE optimization problem of ex22.c with AD for adjoint.\n\n";
4: #include petscda.h
5: #include petscpf.h
6: #include petscmg.h
7: #include petscsnes.h
9: /*
11: Minimize F(w,u) such that G(w,u) = 0
13: L(w,u,lambda) = F(w,u) + lambda^T G(w,u)
15: w - design variables (what we change to get an optimal solution)
16: u - state variables (i.e. the PDE solution)
17: lambda - the Lagrange multipliers
19: U = (w u lambda)
21: fu, fw, flambda contain the gradient of L(w,u,lambda)
23: FU = (fw fu flambda)
25: In this example the PDE is
26: Uxx - u^2 = 2,
27: u(0) = w(0), thus this is the free parameter
28: u(1) = 0
29: the function we wish to minimize is
30: \integral u^{2}
32: The exact solution for u is given by u(x) = x*x - 1.25*x + .25
34: Use the usual centered finite differences.
36: Note we treat the problem as non-linear though it happens to be linear
38: The lambda and u are NOT interlaced.
40: We optionally provide a preconditioner on each level from the operator
42: (1 0 0)
43: (0 J 0)
44: (0 0 J')
46:
47: */
53: typedef struct {
54: Mat J; /* Jacobian of PDE system */
55: KSP ksp; /* Solver for that Jacobian */
56: } AppCtx;
60: PetscErrorCode myPCApply(DMMG dmmg,Vec x,Vec y)
61: {
62: Vec xu,xlambda,yu,ylambda;
63: PetscScalar *xw,*yw;
65: VecPack packer = (VecPack)dmmg->dm;
66: AppCtx *appctx = (AppCtx*)dmmg->user;
69: VecPackGetAccess(packer,x,&xw,&xu,&xlambda);
70: VecPackGetAccess(packer,y,&yw,&yu,&ylambda);
71: if (yw && xw) {
72: yw[0] = xw[0];
73: }
74: KSPSolve(appctx->ksp,xu,yu);
76: KSPSolveTranspose(appctx->ksp,xlambda,ylambda);
77: /* VecCopy(xu,yu);
78: VecCopy(xlambda,ylambda); */
79: VecPackRestoreAccess(packer,x,&xw,&xu,&xlambda);
80: VecPackRestoreAccess(packer,y,&yw,&yu,&ylambda);
81: return(0);
82: }
86: PetscErrorCode myPCView(DMMG dmmg,PetscViewer v)
87: {
89: AppCtx *appctx = (AppCtx*)dmmg->user;
92: KSPView(appctx->ksp,v);
93: return(0);
94: }
98: int main(int argc,char **argv)
99: {
101: PetscInt nlevels,i,j;
102: DA da;
103: DMMG *dmmg;
104: VecPack packer;
105: AppCtx *appctx;
106: ISColoring iscoloring;
107: PetscTruth bdp;
109: PetscInitialize(&argc,&argv,PETSC_NULL,help);
111: /* Hardwire several options; can be changed at command line */
112: PetscOptionsSetValue("-dmmg_grid_sequence",PETSC_NULL);
113: PetscOptionsSetValue("-ksp_type","fgmres");
114: PetscOptionsSetValue("-ksp_max_it","5");
115: PetscOptionsSetValue("-pc_mg_type","full");
116: PetscOptionsSetValue("-mg_coarse_ksp_type","gmres");
117: PetscOptionsSetValue("-mg_levels_ksp_type","gmres");
118: PetscOptionsSetValue("-mg_coarse_ksp_max_it","6");
119: PetscOptionsSetValue("-mg_levels_ksp_max_it","3");
120: PetscOptionsSetValue("-snes_mf_type","wp");
121: PetscOptionsSetValue("-snes_mf_compute_norma","no");
122: PetscOptionsSetValue("-snes_mf_compute_normu","no");
123: PetscOptionsSetValue("-snes_ls","basic");
124: PetscOptionsSetValue("-dmmg_jacobian_mf_fd",0);
125: /* PetscOptionsSetValue("-snes_ls","basicnonorms"); */
126: PetscOptionsInsert(&argc,&argv,PETSC_NULL);
128: /* create VecPack object to manage composite vector */
129: VecPackCreate(PETSC_COMM_WORLD,&packer);
130: VecPackAddArray(packer,1);
131: DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,-5,1,1,PETSC_NULL,&da);
132: VecPackAddDA(packer,da);
133: VecPackAddDA(packer,da);
134: DADestroy(da);
136: /* create nonlinear multi-level solver */
137: DMMGCreate(PETSC_COMM_WORLD,2,PETSC_NULL,&dmmg);
138: DMMGSetDM(dmmg,(DM)packer);
139: VecPackDestroy(packer);
141: /* Create Jacobian of PDE function for each level */
142: nlevels = DMMGGetLevels(dmmg);
143: for (i=0; i<nlevels; i++) {
144: packer = (VecPack)dmmg[i]->dm;
145: VecPackGetEntries(packer,PETSC_NULL,&da,PETSC_NULL);
146: PetscNew(AppCtx,&appctx);
147: DAGetColoring(da,IS_COLORING_GHOSTED,&iscoloring);
148: DAGetMatrix(da,MATAIJ,&appctx->J);
149: MatSetColoring(appctx->J,iscoloring);
150: ISColoringDestroy(iscoloring);
151: DASetLocalFunction(da,(DALocalFunction1)PDEFormFunctionLocal);
152: DASetLocalAdicFunction(da,ad_PDEFormFunctionLocal);
153: dmmg[i]->user = (void*)appctx;
154: }
156: DMMGSetSNES(dmmg,FormFunction,PETSC_NULL);
158: PetscOptionsHasName(PETSC_NULL,"-bdp",&bdp);
159: if (bdp) {
160: for (i=0; i<nlevels; i++) {
161: KSP ksp;
162: PC pc,mpc;
164: appctx = (AppCtx*) dmmg[i]->user;
165: KSPCreate(PETSC_COMM_WORLD,&appctx->ksp);
166: KSPSetOptionsPrefix(appctx->ksp,"bdp_");
167: KSPSetFromOptions(appctx->ksp);
169: SNESGetKSP(dmmg[i]->snes,&ksp);
170: KSPGetPC(ksp,&pc);
171: for (j=0; j<=i; j++) {
172: MGGetSmoother(pc,j,&ksp);
173: KSPGetPC(ksp,&mpc);
174: PCSetType(mpc,PCSHELL);
175: PCShellSetApply(mpc,(PetscErrorCode (*)(void*,Vec,Vec))myPCApply,dmmg[j]);
176: PCShellSetView(mpc,(PetscErrorCode (*)(void*,PetscViewer))myPCView);
177: }
178: }
179: }
181: DMMGSolve(dmmg);
183: /* VecView(DMMGGetx(dmmg),PETSC_VIEWER_SOCKET_WORLD); */
184: for (i=0; i<nlevels; i++) {
185: appctx = (AppCtx*)dmmg[i]->user;
186: MatDestroy(appctx->J);
187: if (appctx->ksp) {KSPDestroy(appctx->ksp);}
188: PetscFree(appctx);
189: }
190: DMMGDestroy(dmmg);
192: PetscFinalize();
193: return 0;
194: }
195:
196: /*
197: Enforces the PDE on the grid
198: This local function acts on the ghosted version of U (accessed via DAGetLocalVector())
199: BUT the global, nonghosted version of FU
201: Process adiC(36): PDEFormFunctionLocal
202: */
205: PetscErrorCode PDEFormFunctionLocal(DALocalInfo *info,PetscScalar *u,PetscScalar *fu,PassiveScalar *w)
206: {
207: PetscInt xs = info->xs,xm = info->xm,i,mx = info->mx;
208: PetscScalar d,h;
210: d = mx-1.0;
211: h = 1.0/d;
213: for (i=xs; i<xs+xm; i++) {
214: if (i == 0) fu[i] = 2.0*d*(u[i] - w[0]) + h*u[i]*u[i];
215: else if (i == mx-1) fu[i] = 2.0*d*u[i] + h*u[i]*u[i];
216: else fu[i] = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h) + h*u[i]*u[i];
217: }
219: PetscLogFlops(9*mx);
220: return 0;
221: }
223: /*
224: Evaluates FU = Gradiant(L(w,u,lambda))
226: This is the function that is usually passed to the SNESSetJacobian() or DMMGSetSNES() and
227: defines the nonlinear set of equations that are to be solved.
229: This local function acts on the ghosted version of U (accessed via VecPackGetLocalVectors() and
230: VecPackScatter()) BUT the global, nonghosted version of FU (via VecPackAccess()).
232: This function uses PDEFormFunction() to enforce the PDE constraint equations and its adjoint
233: for the Lagrange multiplier equations
235: */
238: PetscErrorCode FormFunction(SNES snes,Vec U,Vec FU,void* dummy)
239: {
240: DMMG dmmg = (DMMG)dummy;
242: PetscInt xs,xm,i,N,nredundant;
243: PetscScalar *u,*w,*fw,*fu,*lambda,*flambda,d,h,h2;
244: Vec vu,vlambda,vfu,vflambda,vglambda;
245: DA da;
246: VecPack packer = (VecPack)dmmg->dm;
247: AppCtx *appctx = (AppCtx*)dmmg->user;
248: PetscTruth skipadic;
251: PetscOptionsHasName(0,"-skipadic",&skipadic);
253: VecPackGetEntries(packer,&nredundant,&da,PETSC_IGNORE);
254: DAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
255: DAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0);
256: d = (N-1.0);
257: h = 1.0/d;
258: h2 = 2.0*h;
260: VecPackGetLocalVectors(packer,&w,&vu,&vlambda);
261: VecPackScatter(packer,U,w,vu,vlambda);
262: VecPackGetAccess(packer,FU,&fw,&vfu,&vflambda);
263: VecPackGetAccess(packer,U,0,0,&vglambda);
265: /* G() */
266: DAFormFunction1(da,vu,vfu,w);
267: if (!skipadic) {
268: /* lambda^T G_u() */
269: DAComputeJacobian1WithAdic(da,vu,appctx->J,w);
270: if (appctx->ksp) {
271: KSPSetOperators(appctx->ksp,appctx->J,appctx->J,SAME_NONZERO_PATTERN);
272: }
273: MatMultTranspose(appctx->J,vglambda,vflambda);
274: }
276: DAVecGetArray(da,vu,&u);
277: DAVecGetArray(da,vfu,&fu);
278: DAVecGetArray(da,vlambda,&lambda);
279: DAVecGetArray(da,vflambda,&flambda);
281: /* L_w */
282: if (xs == 0) { /* only first processor computes this */
283: fw[0] = -2.*d*lambda[0];
284: }
286: /* lambda^T G_u() */
287: if (skipadic) {
288: for (i=xs; i<xs+xm; i++) {
289: if (i == 0) flambda[0] = 2.*d*lambda[0] - d*lambda[1] + h2*lambda[0]*u[0];
290: else if (i == 1) flambda[1] = 2.*d*lambda[1] - d*lambda[2] + h2*lambda[1]*u[1];
291: else if (i == N-1) flambda[N-1] = 2.*d*lambda[N-1] - d*lambda[N-2] + h2*lambda[N-1]*u[N-1];
292: else if (i == N-2) flambda[N-2] = 2.*d*lambda[N-2] - d*lambda[N-3] + h2*lambda[N-2]*u[N-2];
293: else flambda[i] = - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]) + h2*lambda[i]*u[i];
294: }
295: }
297: /* F_u */
298: for (i=xs; i<xs+xm; i++) {
299: if (i == 0) flambda[0] += h*u[0];
300: else if (i == 1) flambda[1] += h2*u[1];
301: else if (i == N-1) flambda[N-1] += h*u[N-1];
302: else if (i == N-2) flambda[N-2] += h2*u[N-2];
303: else flambda[i] += h2*u[i];
304: }
306: DAVecRestoreArray(da,vu,&u);
307: DAVecRestoreArray(da,vfu,&fu);
308: DAVecRestoreArray(da,vlambda,&lambda);
309: DAVecRestoreArray(da,vflambda,&flambda);
311: VecPackRestoreLocalVectors(packer,&w,&vu,&vlambda);
312: VecPackRestoreAccess(packer,FU,&fw,&vfu,&vflambda);
313: VecPackRestoreAccess(packer,U,0,0,&vglambda);
315: PetscLogFlops(9*N);
316: return(0);
317: }