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