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
 8:  #include petscdmmg.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: */



 54: typedef struct {
 55:   Mat        J;           /* Jacobian of PDE system */
 56:   KSP       ksp;        /* Solver for that Jacobian */
 57: } AppCtx;

 61: PetscErrorCode myPCApply(DMMG dmmg,Vec x,Vec y)
 62: {
 63:   Vec            xu,xlambda,yu,ylambda;
 64:   PetscScalar    *xw,*yw;
 66:   DMComposite        packer = (DMComposite)dmmg->dm;
 67:   AppCtx         *appctx = (AppCtx*)dmmg->user;

 70:   DMCompositeGetAccess(packer,x,&xw,&xu,&xlambda);
 71:   DMCompositeGetAccess(packer,y,&yw,&yu,&ylambda);
 72:   if (yw && xw) {
 73:     yw[0] = xw[0];
 74:   }
 75:   KSPSolve(appctx->ksp,xu,yu);

 77:   KSPSolveTranspose(appctx->ksp,xlambda,ylambda);
 78:   /*  VecCopy(xu,yu);
 79:       VecCopy(xlambda,ylambda); */
 80:   DMCompositeRestoreAccess(packer,x,&xw,&xu,&xlambda);
 81:   DMCompositeRestoreAccess(packer,y,&yw,&yu,&ylambda);
 82:   return(0);
 83: }

 87: PetscErrorCode myPCView(DMMG dmmg,PetscViewer v)
 88: {
 90:   AppCtx         *appctx = (AppCtx*)dmmg->user;

 93:   KSPView(appctx->ksp,v);
 94:   return(0);
 95: }

 99: int main(int argc,char **argv)
100: {
102:   PetscInt       nlevels,i,j;
103:   DA             da;
104:   DMMG           *dmmg;
105:   DMComposite        packer;
106:   AppCtx         *appctx;
107:   ISColoring     iscoloring;
108:   PetscTruth     bdp;

110:   PetscInitialize(&argc,&argv,PETSC_NULL,help);

112:   /* Hardwire several options; can be changed at command line */
113:   PetscOptionsSetValue("-dmmg_grid_sequence",PETSC_NULL);
114:   PetscOptionsSetValue("-ksp_type","fgmres");
115:   PetscOptionsSetValue("-ksp_max_it","5");
116:   PetscOptionsSetValue("-pc_mg_type","full");
117:   PetscOptionsSetValue("-mg_coarse_ksp_type","gmres");
118:   PetscOptionsSetValue("-mg_levels_ksp_type","gmres");
119:   PetscOptionsSetValue("-mg_coarse_ksp_max_it","6");
120:   PetscOptionsSetValue("-mg_levels_ksp_max_it","3");
121:   PetscOptionsSetValue("-snes_mf_type","wp");
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 DMComposite object to manage composite vector */
129:   DMCompositeCreate(PETSC_COMM_WORLD,&packer);
130:   DMCompositeAddArray(packer,0,1);
131:   DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,-5,1,1,PETSC_NULL,&da);
132:   DMCompositeAddDM(packer,(DM)da);
133:   DMCompositeAddDM(packer,(DM)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:   DMCompositeDestroy(packer);

141:   /* Create Jacobian of PDE function for each level */
142:   nlevels = DMMGGetLevels(dmmg);
143:   for (i=0; i<nlevels; i++) {
144:     packer = (DMComposite)dmmg[i]->dm;
145:     DMCompositeGetEntries(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);
157:   DMMGSetFromOptions(dmmg);

159:   PetscOptionsHasName(PETSC_NULL,"-bdp",&bdp);
160:   if (bdp) {
161:     for (i=0; i<nlevels; i++) {
162:       KSP  ksp;
163:       PC   pc,mpc;

165:       appctx = (AppCtx*) dmmg[i]->user;
166:       KSPCreate(PETSC_COMM_WORLD,&appctx->ksp);
167:       KSPSetOptionsPrefix(appctx->ksp,"bdp_");
168:       KSPSetFromOptions(appctx->ksp);

170:       SNESGetKSP(dmmg[i]->snes,&ksp);
171:       KSPGetPC(ksp,&pc);
172:       for (j=0; j<=i; j++) {
173:         PCMGGetSmoother(pc,j,&ksp);
174:         KSPGetPC(ksp,&mpc);
175:         PCSetType(mpc,PCSHELL);
176:         PCShellSetContext(mpc,dmmg[j]);
177:         PCShellSetApply(mpc,(PetscErrorCode (*)(void*,Vec,Vec))myPCApply);
178:         PCShellSetView(mpc,(PetscErrorCode (*)(void*,PetscViewer))myPCView);
179:       }
180:     }
181:   }

183:   DMMGSolve(dmmg);

185:   for (i=0; i<nlevels; i++) {
186:     appctx = (AppCtx*)dmmg[i]->user;
187:     MatDestroy(appctx->J);
188:     if (appctx->ksp) {KSPDestroy(appctx->ksp);}
189:     PetscFree(appctx);
190:   }
191:   DMMGDestroy(dmmg);

193:   PetscFinalize();
194:   return 0;
195: }
196: 
197: /*
198:      Enforces the PDE on the grid
199:      This local function acts on the ghosted version of U (accessed via DAGetLocalVector())
200:      BUT the global, nonghosted version of FU

202:      Process adiC(36): PDEFormFunctionLocal
203: */
206: PetscErrorCode PDEFormFunctionLocal(DALocalInfo *info,PetscScalar *u,PetscScalar *fu,PassiveScalar *w)
207: {
208:   PetscInt       xs = info->xs,xm = info->xm,i,mx = info->mx;
209:   PetscScalar    d,h;

212:   d    = mx-1.0;
213:   h    = 1.0/d;

215:   for (i=xs; i<xs+xm; i++) {
216:     if      (i == 0)    fu[i]   = 2.0*d*(u[i] - w[0]) + h*u[i]*u[i];
217:     else if (i == mx-1) fu[i]   = 2.0*d*u[i] + h*u[i]*u[i];
218:     else                fu[i]   = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h) + h*u[i]*u[i];
219:   }

221:   PetscLogFlops(9*mx);
222:   return 0;
223: }

225: /*
226:       Evaluates FU = Gradiant(L(w,u,lambda))

228:       This is the function that is usually passed to the SNESSetJacobian() or DMMGSetSNES() and
229:     defines the nonlinear set of equations that are to be solved.

231:      This local function acts on the ghosted version of U (accessed via DMCompositeGetLocalVectors() and
232:    DMCompositeScatter()) BUT the global, nonghosted version of FU (via DMCompositeAccess()).

234:      This function uses PDEFormFunction() to enforce the PDE constraint equations and its adjoint
235:    for the Lagrange multiplier equations

237: */
240: PetscErrorCode FormFunction(SNES snes,Vec U,Vec FU,void* dummy)
241: {
242:   DMMG           dmmg = (DMMG)dummy;
244:   PetscInt       xs,xm,i,N,nredundant;
245:   PetscScalar    *u,*w,*fw,*fu,*lambda,*flambda,d,h,h2;
246:   Vec            vu,vlambda,vfu,vflambda,vglambda;
247:   DA             da;
248:   DMComposite        packer = (DMComposite)dmmg->dm;
249:   PetscTruth     useadic = PETSC_TRUE;
250: #if defined(PETSC_HAVE_ADIC)
251:   AppCtx         *appctx = (AppCtx*)dmmg->user;
252: #endif


256: #if defined(PETSC_HAVE_ADIC)
257:   PetscOptionsHasName(0,"-useadic",&skipadic);
258: #endif

260:   DMCompositeGetEntries(packer,&nredundant,&da,PETSC_IGNORE);
261:   DAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
262:   DAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0);
263:   d    = (N-1.0);
264:   h    = 1.0/d;
265:   h2   = 2.0*h;

267:   DMCompositeGetLocalVectors(packer,&w,&vu,&vlambda);
268:   DMCompositeScatter(packer,U,w,vu,vlambda);
269:   DMCompositeGetAccess(packer,FU,&fw,&vfu,&vflambda);
270:   DMCompositeGetAccess(packer,U,0,0,&vglambda);

272:   /* G() */
273:   DAFormFunction1(da,vu,vfu,w);

275: #if defined(PETSC_HAVE_ADIC)
276:   if (useadic) {
277:     /* lambda^T G_u() */
278:     DAComputeJacobian1WithAdic(da,vu,appctx->J,w);
279:     if (appctx->ksp) {
280:       KSPSetOperators(appctx->ksp,appctx->J,appctx->J,SAME_NONZERO_PATTERN);
281:     }
282:     MatMultTranspose(appctx->J,vglambda,vflambda);
283:   }
284: #endif

286:   DAVecGetArray(da,vu,&u);
287:   DAVecGetArray(da,vfu,&fu);
288:   DAVecGetArray(da,vlambda,&lambda);
289:   DAVecGetArray(da,vflambda,&flambda);

291:   /* L_w */
292:   if (xs == 0) { /* only first processor computes this */
293:     fw[0] = -2.*d*lambda[0];
294:   }

296:   /* lambda^T G_u() */
297:   if (!useadic) {
298:     for (i=xs; i<xs+xm; i++) {
299:       if      (i == 0)   flambda[0]   = 2.*d*lambda[0]   - d*lambda[1] + h2*lambda[0]*u[0];
300:       else if (i == 1)   flambda[1]   = 2.*d*lambda[1]   - d*lambda[2] + h2*lambda[1]*u[1];
301:       else if (i == N-1) flambda[N-1] = 2.*d*lambda[N-1] - d*lambda[N-2] + h2*lambda[N-1]*u[N-1];
302:       else if (i == N-2) flambda[N-2] = 2.*d*lambda[N-2] - d*lambda[N-3] + h2*lambda[N-2]*u[N-2];
303:       else               flambda[i]   = - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]) + h2*lambda[i]*u[i];
304:     }
305:   }

307:   /* F_u */
308:   for (i=xs; i<xs+xm; i++) {
309:     if      (i == 0)   flambda[0]   +=    h*u[0];
310:     else if (i == 1)   flambda[1]   +=    h2*u[1];
311:     else if (i == N-1) flambda[N-1] +=    h*u[N-1];
312:     else if (i == N-2) flambda[N-2] +=    h2*u[N-2];
313:     else               flambda[i]   +=    h2*u[i];
314:   }

316:   DAVecRestoreArray(da,vu,&u);
317:   DAVecRestoreArray(da,vfu,&fu);
318:   DAVecRestoreArray(da,vlambda,&lambda);
319:   DAVecRestoreArray(da,vflambda,&flambda);

321:   DMCompositeRestoreLocalVectors(packer,&w,&vu,&vlambda);
322:   DMCompositeRestoreAccess(packer,FU,&fw,&vfu,&vflambda);
323:   DMCompositeRestoreAccess(packer,U,0,0,&vglambda);

325:   PetscLogFlops(9*N);
326:   return(0);
327: }