Actual source code: ex22.c

  1: /*$Id: ex22.c,v 1.18 2001/04/10 19:37:05 bsmith Exp $*/

  3: static char help[] = "Solves PDE optimization problem.nn";

 5:  #include petscda.h
 6:  #include petscpf.h
 7:  #include petscsnes.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;

 44: extern int FormFunction(SNES,Vec,Vec,void*);
 45: extern int Monitor(SNES,int,PetscReal,void*);


 48: int main(int argc,char **argv)
 49: {
 50:   int     ierr,N = 5,i;
 51:   UserCtx user;
 52:   DA      da;
 53:   DMMG    *dmmg;
 54:   VecPack packer;

 56:   PetscInitialize(&argc,&argv,PETSC_NULL,help);
 57:   PetscOptionsGetInt(PETSC_NULL,"-N",&N,PETSC_NULL);

 59:   /* Hardwire several options; can be changed at command line */
 60:   PetscOptionsSetValue("-dmmg_grid_sequence",PETSC_NULL);
 61:   PetscOptionsSetValue("-ksp_type","fgmres");
 62:   PetscOptionsSetValue("-ksp_max_it","5");
 63:   PetscOptionsSetValue("-pc_mg_type","full");
 64:   PetscOptionsSetValue("-mg_coarse_ksp_type","gmres");
 65:   PetscOptionsSetValue("-mg_levels_ksp_type","gmres");
 66:   PetscOptionsSetValue("-mg_coarse_ksp_max_it","6");
 67:   PetscOptionsSetValue("-mg_levels_ksp_max_it","3");
 68:   PetscOptionsSetValue("-snes_mf_type","wp");
 69:   PetscOptionsSetValue("-snes_mf_compute_norma","no");
 70:   PetscOptionsSetValue("-snes_mf_compute_normu","no");
 71:   PetscOptionsSetValue("-snes_eq_ls","basic");
 72:   /* PetscOptionsSetValue("-snes_eq_ls","basicnonorms"); */
 73:   PetscOptionsInsert(&argc,&argv,PETSC_NULL);
 74: 
 75:   /* Create a global vector that includes a single redundant array and two da arrays */
 76:   VecPackCreate(PETSC_COMM_WORLD,&packer);
 77:   VecPackAddArray(packer,1);
 78:   DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,N,2,1,PETSC_NULL,&da);
 79:   VecPackAddDA(packer,da);

 81:   /* create graphics windows */
 82:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u_lambda - state variables and Lagrange multipliers",-1,-1,-1,-1,&user.u_lambda_viewer);
 83:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu_lambda - derivate w.r.t. state variables and Lagrange multipliers",-1,-1,-1,-1,&user.fu_lambda_viewer);

 85:   /* create nonlinear multi-level solver */
 86:   DMMGCreate(PETSC_COMM_WORLD,2,&user,&dmmg);
 87:   DMMGSetUseMatrixFree(dmmg);
 88:   DMMGSetDM(dmmg,(DM)packer);
 89:   DMMGSetSNES(dmmg,FormFunction,PETSC_NULL);
 90:   for (i=0; i<DMMGGetLevels(dmmg); i++) {
 91:     SNESSetMonitor(dmmg[i]->snes,Monitor,dmmg[i],0);
 92:   }
 93:   DMMGSolve(dmmg);
 94:   DMMGDestroy(dmmg);

 96:   DADestroy(da);
 97:   VecPackDestroy(packer);
 98:   PetscViewerDestroy(user.u_lambda_viewer);
 99:   PetscViewerDestroy(user.fu_lambda_viewer);

101:   PetscFinalize();
102:   return 0;
103: }

105: typedef struct {
106:   Scalar u;
107:   Scalar lambda;
108: } ULambda;
109: 
110: /*
111:       Evaluates FU = Gradiant(L(w,u,lambda))

113:      This local function acts on the ghosted version of U (accessed via VecPackGetLocalVectors() and
114:    VecPackScatter()) BUT the global, nonghosted version of FU (via VecPackGetAccess()).

116: */
117: int FormFunction(SNES snes,Vec U,Vec FU,void* dummy)
118: {
119:   DMMG    dmmg = (DMMG)dummy;
120:   int     ierr,xs,xm,i,N,nredundant;
121:   ULambda *u_lambda,*fu_lambda;
122:   Scalar  d,h,*w,*fw;
123:   Vec     vu_lambda,vfu_lambda;
124:   DA      da;
125:   VecPack packer = (VecPack)dmmg->dm;

128:   VecPackGetEntries(packer,&nredundant,&da);
129:   VecPackGetLocalVectors(packer,&w,&vu_lambda);
130:   VecPackScatter(packer,U,w,vu_lambda);
131:   VecPackGetAccess(packer,FU,&fw,&vfu_lambda);

133:   DAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
134:   DAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0);
135:   DAVecGetArray(da,vu_lambda,(void**)&u_lambda);
136:   DAVecGetArray(da,vfu_lambda,(void**)&fu_lambda);
137:   d    = N-1.0;
138:   h    = 1.0/d;

140: #define u(i)        u_lambda[i][0]
141: #define lambda(i)   u_lambda[i][1]
142: #define fu(i)       fu_lambda[i][1]
143: #define flambda(i)  fu_lambda[i][0]

145:   /* derivative of L() w.r.t. w */
146:   if (xs == 0) { /* only first processor computes this */
147:     fw[0] = -2.0*d*u_lambda[0].lambda;
148:   }

150:   /* derivative of L() w.r.t. u */
151:   for (i=xs; i<xs+xm; i++) {
152:     if      (i == 0)   fu_lambda[0].lambda   =    h*u_lambda[0].u   + 2.*d*u_lambda[0].lambda   - d*u_lambda[1].lambda;
153:     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;
154:     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;
155:     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;
156:     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);
157:   }

159:   /* derivative of L() w.r.t. lambda */
160:   for (i=xs; i<xs+xm; i++) {
161:     if      (i == 0)   fu_lambda[0].u   = 2.0*d*(u_lambda[0].u - w[0]);
162:     else if (i == N-1) fu_lambda[N-1].u = 2.0*d*u_lambda[N-1].u;
163:     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);
164:   }

166:   DAVecRestoreArray(da,vu_lambda,(void**)&u_lambda);
167:   DAVecRestoreArray(da,vfu_lambda,(void**)&fu_lambda);
168:   VecPackRestoreLocalVectors(packer,&w,&vu_lambda);
169:   VecPackRestoreAccess(packer,FU,&fw,&vfu_lambda);
170:   PetscLogFlops(13*N);
171:   return(0);
172: }

174: /* 
175:     Computes the exact solution
176: */
177: int u_solution(void *dummy,int n,Scalar *x,Scalar *u)
178: {
179:   int i;
181:   for (i=0; i<n; i++) {
182:     u[2*i] = x[i]*x[i] - 1.25*x[i] + .25;
183:   }
184:   return(0);
185: }

187: int ExactSolution(VecPack packer,Vec U)
188: {
189:   PF      pf;
190:   Vec     x;
191:   Vec     u_global;
192:   Scalar  *w;
193:   DA      da;
194:   int     m,ierr;

197:   VecPackGetEntries(packer,&m,&da);

199:   PFCreate(PETSC_COMM_WORLD,1,1,&pf);
200:   PFSetType(pf,PFQUICK,(void*)u_solution);
201:   DAGetCoordinates(da,&x);
202:   if (!x) {
203:     DASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);
204:     DAGetCoordinates(da,&x);
205:   }
206:   VecPackGetAccess(packer,U,&w,&u_global,0);
207:   if (w) w[0] = .25;
208:   PFApplyVec(pf,x,u_global);
209:   PFDestroy(pf);
210:   VecPackRestoreAccess(packer,U,&w,&u_global,0);
211:   return(0);
212: }


215: int Monitor(SNES snes,int its,PetscReal rnorm,void *dummy)
216: {
217:   DMMG      dmmg = (DMMG)dummy;
218:   UserCtx   *user = (UserCtx*)dmmg->user;
219:   int       ierr,m,N;
220:   Scalar    mone = -1.0,*w,*dw;
221:   Vec       u_lambda,U,F,Uexact;
222:   VecPack   packer = (VecPack)dmmg->dm;
223:   PetscReal norm;
224:   DA        da;

227:   SNESGetSolution(snes,&U);
228:   VecPackGetAccess(packer,U,&w,&u_lambda);
229:   VecView(u_lambda,user->u_lambda_viewer);
230:   VecPackRestoreAccess(packer,U,&w,&u_lambda);

232:   SNESGetFunction(snes,&F,0,0);
233:   VecPackGetAccess(packer,F,&w,&u_lambda);
234:   /* VecView(u_lambda,user->fu_lambda_viewer); */
235:   VecPackRestoreAccess(packer,U,&w,&u_lambda);

237:   VecPackGetEntries(packer,&m,&da);
238:   DAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0);
239:   VecDuplicate(U,&Uexact);
240:   ExactSolution(packer,Uexact);
241:   VecAXPY(&mone,U,Uexact);
242:   VecPackGetAccess(packer,Uexact,&dw,&u_lambda);
243:   VecStrideNorm(u_lambda,0,NORM_2,&norm);
244:   norm = norm/sqrt(N-1.);
245:   if (dw) PetscPrintf(dmmg->comm,"Norm of error %g Error at x = 0 %gn",norm,dw[0]);
246:   VecView(u_lambda,user->fu_lambda_viewer);
247:   VecPackRestoreAccess(packer,Uexact,&dw,&u_lambda);
248:   VecDestroy(Uexact);
249:   return(0);
250: }