Actual source code: ex24.c
2: static char help[] = "Solves PDE optimization problem of ex22.c with AD for adjoint.\n\n";
4: #include <petscdmda.h>
5: #include <petscdmcomposite.h>
6: #include <petscpf.h>
7: #include <petscpcmg.h>
8: #include <petscsnes.h>
9: #include <petscdmmg.h>
11: /*
13: Minimize F(w,u) such that G(w,u) = 0
15: L(w,u,lambda) = F(w,u) + lambda^T G(w,u)
17: w - design variables (what we change to get an optimal solution)
18: u - state variables (i.e. the PDE solution)
19: lambda - the Lagrange multipliers
21: U = (w u lambda)
23: fu, fw, flambda contain the gradient of L(w,u,lambda)
25: FU = (fw fu flambda)
27: In this example the PDE is
28: Uxx - u^2 = 2,
29: u(0) = w(0), thus this is the free parameter
30: u(1) = 0
31: the function we wish to minimize is
32: \integral u^{2}
34: The exact solution for u is given by u(x) = x*x - 1.25*x + .25
36: Use the usual centered finite differences.
38: Note we treat the problem as non-linear though it happens to be linear
40: The lambda and u are NOT interlaced.
42: We optionally provide a preconditioner on each level from the operator
44: (1 0 0)
45: (0 J 0)
46: (0 0 J')
48:
49: */
55: typedef struct {
56: Mat J; /* Jacobian of PDE system */
57: KSP ksp; /* Solver for that Jacobian */
58: } AppCtx;
62: PetscErrorCode myPCApply(PC pc,Vec x,Vec y)
63: {
64: Vec xu,xlambda,yu,ylambda;
65: PetscScalar *xw,*yw;
67: DMMG dmmg;
68: DM packer;
69: AppCtx *appctx;
72: PCShellGetContext(pc,(void**)&dmmg);
73: packer = dmmg->dm;
74: appctx = (AppCtx*)dmmg->user;
75: DMCompositeGetAccess(packer,x,&xw,&xu,&xlambda);
76: DMCompositeGetAccess(packer,y,&yw,&yu,&ylambda);
77: if (yw && xw) {
78: yw[0] = xw[0];
79: }
80: KSPSolve(appctx->ksp,xu,yu);
82: KSPSolveTranspose(appctx->ksp,xlambda,ylambda);
83: /* VecCopy(xu,yu);
84: VecCopy(xlambda,ylambda); */
85: DMCompositeRestoreAccess(packer,x,&xw,&xu,&xlambda);
86: DMCompositeRestoreAccess(packer,y,&yw,&yu,&ylambda);
87: return(0);
88: }
92: PetscErrorCode myPCView(PC pc,PetscViewer v)
93: {
95: DMMG dmmg;
96: AppCtx *appctx;
99: PCShellGetContext(pc,(void**)&dmmg);
100: appctx = (AppCtx*)dmmg->user;
101: KSPView(appctx->ksp,v);
102: return(0);
103: }
107: int main(int argc,char **argv)
108: {
110: PetscInt nlevels,i,j;
111: DM da;
112: DMMG *dmmg;
113: DM packer;
114: AppCtx *appctx;
115: ISColoring iscoloring;
116: PetscBool bdp;
118: PetscInitialize(&argc,&argv,PETSC_NULL,help);
120: /* Hardwire several options; can be changed at command line */
121: PetscOptionsSetValue("-dmmg_grid_sequence",PETSC_NULL);
122: PetscOptionsSetValue("-ksp_type","fgmres");
123: PetscOptionsSetValue("-ksp_max_it","5");
124: PetscOptionsSetValue("-pc_mg_type","full");
125: PetscOptionsSetValue("-mg_coarse_ksp_type","gmres");
126: PetscOptionsSetValue("-mg_levels_ksp_type","gmres");
127: PetscOptionsSetValue("-mg_coarse_ksp_max_it","6");
128: PetscOptionsSetValue("-mg_levels_ksp_max_it","3");
129: PetscOptionsSetValue("-mat_mffd_type","wp");
130: PetscOptionsSetValue("-mat_mffd_compute_normu","no");
131: PetscOptionsSetValue("-snes_ls","basic");
132: PetscOptionsSetValue("-dmmg_jacobian_mf_fd",0);
133: /* PetscOptionsSetValue("-snes_ls","basicnonorms"); */
134: PetscOptionsInsert(&argc,&argv,PETSC_NULL);
136: /* create DMComposite object to manage composite vector */
137: DMCompositeCreate(PETSC_COMM_WORLD,&packer);
138: DMCompositeAddArray(packer,0,1);
139: DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-5,1,1,PETSC_NULL,&da);
140: DMCompositeAddDM(packer,(DM)da);
141: DMCompositeAddDM(packer,(DM)da);
142: DMDestroy(&da);
144: /* create nonlinear multi-level solver */
145: DMMGCreate(PETSC_COMM_WORLD,2,PETSC_NULL,&dmmg);
146: DMMGSetDM(dmmg,(DM)packer);
147: DMDestroy(&packer);
149: /* Create Jacobian of PDE function for each level */
150: nlevels = DMMGGetLevels(dmmg);
151: for (i=0; i<nlevels; i++) {
152: packer = dmmg[i]->dm;
153: DMCompositeGetEntries(packer,PETSC_NULL,&da,PETSC_NULL);
154: PetscNew(AppCtx,&appctx);
155: DMGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);
156: DMGetMatrix(da,MATAIJ,&appctx->J);
157: MatSetColoring(appctx->J,iscoloring);
158: ISColoringDestroy(&iscoloring);
159: DMDASetLocalFunction(da,(DMDALocalFunction1)PDEFormFunctionLocal);
160: DMDASetLocalAdicFunction(da,ad_PDEFormFunctionLocal);
161: dmmg[i]->user = (void*)appctx;
162: }
164: DMMGSetSNES(dmmg,FormFunction,PETSC_NULL);
165: DMMGSetFromOptions(dmmg);
167: PetscOptionsHasName(PETSC_NULL,"-bdp",&bdp);
168: if (bdp) {
169: for (i=0; i<nlevels; i++) {
170: KSP ksp;
171: PC pc,mpc;
173: appctx = (AppCtx*) dmmg[i]->user;
174: KSPCreate(PETSC_COMM_WORLD,&appctx->ksp);
175: KSPSetOptionsPrefix(appctx->ksp,"bdp_");
176: KSPSetFromOptions(appctx->ksp);
178: SNESGetKSP(dmmg[i]->snes,&ksp);
179: KSPGetPC(ksp,&pc);
180: for (j=0; j<=i; j++) {
181: PCMGGetSmoother(pc,j,&ksp);
182: KSPGetPC(ksp,&mpc);
183: PCSetType(mpc,PCSHELL);
184: PCShellSetContext(mpc,dmmg[j]);
185: PCShellSetApply(mpc,myPCApply);
186: PCShellSetView(mpc,myPCView);
187: }
188: }
189: }
191: DMMGSolve(dmmg);
193: for (i=0; i<nlevels; i++) {
194: appctx = (AppCtx*)dmmg[i]->user;
195: MatDestroy(&appctx->J);
196: if (appctx->ksp) {KSPDestroy(&appctx->ksp);}
197: PetscFree(appctx);
198: }
199: DMMGDestroy(dmmg);
201: PetscFinalize();
202: return 0;
203: }
204:
205: /*
206: Enforces the PDE on the grid
207: This local function acts on the ghosted version of U (accessed via DMGetLocalVector())
208: BUT the global, nonghosted version of FU
210: Process adiC(36): PDEFormFunctionLocal
211: */
214: PetscErrorCode PDEFormFunctionLocal(DMDALocalInfo *info,PetscScalar *u,PetscScalar *fu,PassiveScalar *w)
215: {
216: PetscInt xs = info->xs,xm = info->xm,i,mx = info->mx;
217: PetscScalar d,h;
220: d = mx-1.0;
221: h = 1.0/d;
223: for (i=xs; i<xs+xm; i++) {
224: if (i == 0) fu[i] = 2.0*d*(u[i] - w[0]) + h*u[i]*u[i];
225: else if (i == mx-1) fu[i] = 2.0*d*u[i] + h*u[i]*u[i];
226: else fu[i] = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h) + h*u[i]*u[i];
227: }
229: PetscLogFlops(9.0*mx);
230: return 0;
231: }
233: /*
234: Evaluates FU = Gradiant(L(w,u,lambda))
236: This is the function that is usually passed to the SNESSetJacobian() or DMMGSetSNES() and
237: defines the nonlinear set of equations that are to be solved.
239: This local function acts on the ghosted version of U (accessed via DMCompositeGetLocalVectors() and
240: DMCompositeScatter()) BUT the global, nonghosted version of FU (via DMCompositeAccess()).
242: This function uses PDEFormFunction() to enforce the PDE constraint equations and its adjoint
243: for the Lagrange multiplier equations
245: */
248: PetscErrorCode FormFunction(SNES snes,Vec U,Vec FU,void* dummy)
249: {
250: DMMG dmmg = (DMMG)dummy;
252: PetscInt xs,xm,i,N,nredundant;
253: PetscScalar *u,*w,*fw,*fu,*lambda,*flambda,d,h,h2;
254: Vec vu,vlambda,vfu,vflambda,vglambda;
255: DM da;
256: DM packer = (DM)dmmg->dm;
257: PetscBool useadic = PETSC_TRUE;
258: #if defined(PETSC_HAVE_ADIC)
259: AppCtx *appctx = (AppCtx*)dmmg->user;
260: #endif
264: #if defined(PETSC_HAVE_ADIC)
265: PetscOptionsHasName(0,"-useadic",&skipadic);
266: #endif
268: DMCompositeGetEntries(packer,&nredundant,&da,PETSC_IGNORE);
269: DMDAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
270: DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0);
271: d = (N-1.0);
272: h = 1.0/d;
273: h2 = 2.0*h;
275: DMCompositeGetLocalVectors(packer,&w,&vu,&vlambda);
276: DMCompositeScatter(packer,U,w,vu,vlambda);
277: DMCompositeGetAccess(packer,FU,&fw,&vfu,&vflambda);
278: DMCompositeGetAccess(packer,U,0,0,&vglambda);
280: /* G() */
281: DMDAFormFunction1(da,vu,vfu,w);
283: #if defined(PETSC_HAVE_ADIC)
284: if (useadic) {
285: /* lambda^T G_u() */
286: DMDAComputeJacobian1WithAdic(da,vu,appctx->J,w);
287: if (appctx->ksp) {
288: KSPSetOperators(appctx->ksp,appctx->J,appctx->J,SAME_NONZERO_PATTERN);
289: }
290: MatMultTranspose(appctx->J,vglambda,vflambda);
291: }
292: #endif
294: DMDAVecGetArray(da,vu,&u);
295: DMDAVecGetArray(da,vfu,&fu);
296: DMDAVecGetArray(da,vlambda,&lambda);
297: DMDAVecGetArray(da,vflambda,&flambda);
299: /* L_w */
300: if (xs == 0) { /* only first processor computes this */
301: fw[0] = -2.*d*lambda[0];
302: }
304: /* lambda^T G_u() */
305: if (!useadic) {
306: for (i=xs; i<xs+xm; i++) {
307: if (i == 0) flambda[0] = 2.*d*lambda[0] - d*lambda[1] + h2*lambda[0]*u[0];
308: else if (i == 1) flambda[1] = 2.*d*lambda[1] - d*lambda[2] + h2*lambda[1]*u[1];
309: else if (i == N-1) flambda[N-1] = 2.*d*lambda[N-1] - d*lambda[N-2] + h2*lambda[N-1]*u[N-1];
310: else if (i == N-2) flambda[N-2] = 2.*d*lambda[N-2] - d*lambda[N-3] + h2*lambda[N-2]*u[N-2];
311: else flambda[i] = - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]) + h2*lambda[i]*u[i];
312: }
313: }
315: /* F_u */
316: for (i=xs; i<xs+xm; i++) {
317: if (i == 0) flambda[0] += h*u[0];
318: else if (i == 1) flambda[1] += h2*u[1];
319: else if (i == N-1) flambda[N-1] += h*u[N-1];
320: else if (i == N-2) flambda[N-2] += h2*u[N-2];
321: else flambda[i] += h2*u[i];
322: }
324: DMDAVecRestoreArray(da,vu,&u);
325: DMDAVecRestoreArray(da,vfu,&fu);
326: DMDAVecRestoreArray(da,vlambda,&lambda);
327: DMDAVecRestoreArray(da,vflambda,&flambda);
329: DMCompositeRestoreLocalVectors(packer,&w,&vu,&vlambda);
330: DMCompositeRestoreAccess(packer,FU,&fw,&vfu,&vflambda);
331: DMCompositeRestoreAccess(packer,U,0,0,&vglambda);
333: PetscLogFlops(9.0*N);
334: return(0);
335: }