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: #undef __FUNCT__
 50: int main(int argc,char **argv)
 51: {
 52:   int     ierr;
 53:   UserCtx user;
 54:   DA      da;
 55:   DMMG    *dmmg;
 56:   VecPack packer;

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

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

 77:   /* Create a global vector that includes a single redundant array and two da arrays */
 78:   VecPackCreate(PETSC_COMM_WORLD,&packer);
 79:   VecPackAddArray(packer,1);
 80:   DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,-5,2,1,PETSC_NULL,&da);
 81:   VecPackAddDA(packer,da);

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

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

 98:   DADestroy(da);
 99:   VecPackDestroy(packer);
100:   PetscViewerDestroy(user.u_lambda_viewer);
101:   PetscViewerDestroy(user.fu_lambda_viewer);

103:   PetscFinalize();
104:   return 0;
105: }

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

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

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