Actual source code: ex22.c

  1: /*$Id: ex22.c,v 1.23 2001/08/07 21:31:17 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;
 51:   UserCtx user;
 52:   DA      da;
 53:   DMMG    *dmmg;
 54:   VecPack packer;

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

 58:   /* Hardwire several options; can be changed at command line */
 59:   PetscOptionsSetValue("-dmmg_grid_sequence",PETSC_NULL);
 60:   PetscOptionsSetValue("-ksp_type","fgmres");
 61:   PetscOptionsSetValue("-ksp_max_it","5");
 62:   PetscOptionsSetValue("-pc_mg_type","full");
 63:   PetscOptionsSetValue("-mg_coarse_ksp_type","gmres");
 64:   PetscOptionsSetValue("-mg_levels_ksp_type","gmres");
 65:   PetscOptionsSetValue("-mg_coarse_ksp_max_it","6");
 66:   PetscOptionsSetValue("-mg_levels_ksp_max_it","3");
 67:   PetscOptionsSetValue("-snes_mf_type","wp");
 68:   PetscOptionsSetValue("-snes_mf_compute_norma","no");
 69:   PetscOptionsSetValue("-snes_mf_compute_normu","no");
 70:   PetscOptionsSetValue("-snes_ls","basic");
 71:   PetscOptionsSetValue("-dmmg_jacobian_mf_fd",0);
 72:   /* PetscOptionsSetValue("-snes_ls","basicnonorms"); */
 73:   PetscOptionsInsert(&argc,&argv,PETSC_NULL);

 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,-5,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:   DMMGSetDM(dmmg,(DM)packer);
 88:   DMMGSetSNES(dmmg,FormFunction,PETSC_NULL);
 89:   /*
 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:   PetscScalar u;
107:   PetscScalar 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:   PetscScalar  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:   /* derivative of L() w.r.t. w */
141:   if (xs == 0) { /* only first processor computes this */
142:     fw[0] = -2.0*d*u_lambda[0].lambda;
143:   }

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

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

161:   DAVecRestoreArray(da,vu_lambda,(void**)&u_lambda);
162:   DAVecRestoreArray(da,vfu_lambda,(void**)&fu_lambda);
163:   VecPackRestoreLocalVectors(packer,&w,&vu_lambda);
164:   VecPackRestoreAccess(packer,FU,&fw,&vfu_lambda);
165:   PetscLogFlops(13*N);
166:   return(0);
167: }

169: /* 
170:     Computes the exact solution
171: */
172: int u_solution(void *dummy,int n,PetscScalar *x,PetscScalar *u)
173: {
174:   int i;
176:   for (i=0; i<n; i++) {
177:     u[2*i] = x[i]*x[i] - 1.25*x[i] + .25;
178:   }
179:   return(0);
180: }

182: int ExactSolution(VecPack packer,Vec U)
183: {
184:   PF      pf;
185:   Vec     x;
186:   Vec     u_global;
187:   PetscScalar  *w;
188:   DA      da;
189:   int     m,ierr;

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

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


210: int Monitor(SNES snes,int its,PetscReal rnorm,void *dummy)
211: {
212:   DMMG         dmmg = (DMMG)dummy;
213:   UserCtx      *user = (UserCtx*)dmmg->user;
214:   int          ierr,m,N;
215:   PetscScalar  mone = -1.0,*w,*dw;
216:   Vec          u_lambda,U,F,Uexact;
217:   VecPack      packer = (VecPack)dmmg->dm;
218:   PetscReal    norm;
219:   DA           da;

222:   SNESGetSolution(snes,&U);
223:   VecPackGetAccess(packer,U,&w,&u_lambda);
224:   VecView(u_lambda,user->u_lambda_viewer);
225:   VecPackRestoreAccess(packer,U,&w,&u_lambda);

227:   SNESGetFunction(snes,&F,0,0);
228:   VecPackGetAccess(packer,F,&w,&u_lambda);
229:   /* VecView(u_lambda,user->fu_lambda_viewer); */
230:   VecPackRestoreAccess(packer,U,&w,&u_lambda);

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