Actual source code: ex22.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: static const char help[] = "Solves PDE optimization problem using full-space method, interlaces state and adjoint variables.\n\n";

  4: #include <petscdm.h>
  5: #include <petscdmda.h>
  6: #include <petscdmredundant.h>
  7: #include <petscdmcomposite.h>
  8: #include <petscpf.h>
  9: #include <petscsnes.h>
 10: #include <petsc/private/dmimpl.h>

 12: /*

 14:        w - design variables (what we change to get an optimal solution)
 15:        u - state variables (i.e. the PDE solution)
 16:        lambda - the Lagrange multipliers

 18:             U = (w [u_0 lambda_0 u_1 lambda_1 .....])

 20:        fu, fw, flambda contain the gradient of L(w,u,lambda)

 22:             FU = (fw [fu_0 flambda_0 .....])

 24:        In this example the PDE is
 25:                              Uxx = 2,
 26:                             u(0) = w(0), thus this is the free parameter
 27:                             u(1) = 0
 28:        the function we wish to minimize is
 29:                             \integral u^{2}

 31:        The exact solution for u is given by u(x) = x*x - 1.25*x + .25

 33:        Use the usual centered finite differences.

 35:        Note we treat the problem as non-linear though it happens to be linear

 37:        See ex21.c for the same code, but that does NOT interlaces the u and the lambda

 39:        The vectors u_lambda and fu_lambda contain the u and the lambda interlaced
 40: */

 42: typedef struct {
 43:   PetscViewer u_lambda_viewer;
 44:   PetscViewer fu_lambda_viewer;
 45: } UserCtx;

 47: extern PetscErrorCode ComputeFunction(SNES,Vec,Vec,void*);
 48: extern PetscErrorCode ComputeJacobian_MF(SNES,Vec,Mat,Mat,void*);
 49: extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);

 51: /*
 52:     Uses full multigrid preconditioner with GMRES (with no preconditioner inside the GMRES) as the
 53:   smoother on all levels. This is because (1) in the matrix free case no matrix entries are
 54:   available for doing Jacobi or SOR preconditioning and (2) the explicit matrix case the diagonal
 55:   entry for the control variable is zero which means default SOR will not work.

 57: */
 58: char common_options[] = "-ksp_type fgmres\
 59:                          -snes_grid_sequence 4 \
 60:                          -pc_type mg\
 61:                          -mg_levels_pc_type none \
 62:                          -mg_coarse_pc_type none \
 63:                          -pc_mg_type full \
 64:                          -mg_coarse_ksp_type gmres \
 65:                          -mg_levels_ksp_type gmres \
 66:                          -mg_coarse_ksp_max_it 6 \
 67:                          -mg_levels_ksp_max_it 3";

 69: char matrix_free_options[] = "-mat_mffd_compute_normu no \
 70:                               -mat_mffd_type wp";

 72: extern PetscErrorCode DMCreateMatrix_MF(DM,Mat*);

 76: int main(int argc,char **argv)
 77: {
 79:   UserCtx        user;
 80:   DM             red,da;
 81:   SNES           snes;
 82:   DM             packer;
 83:   PetscBool      use_monitor = PETSC_FALSE;

 85:   PetscInitialize(&argc,&argv,NULL,help);
 86:   PetscOptionsSetFromOptions();

 88:   /* Hardwire several options; can be changed at command line */
 89:   PetscOptionsInsertString(common_options);
 90:   PetscOptionsInsertString(matrix_free_options);
 91:   PetscOptionsInsert(&argc,&argv,NULL);
 92:   PetscOptionsGetBool(NULL,"-use_monitor",&use_monitor,PETSC_IGNORE);

 94:   /* Create a global vector that includes a single redundant array and two da arrays */
 95:   DMCompositeCreate(PETSC_COMM_WORLD,&packer);
 96:   DMRedundantCreate(PETSC_COMM_WORLD,0,1,&red);
 97:   DMSetOptionsPrefix(red,"red_");
 98:   DMCompositeAddDM(packer,red);
 99:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,-5,2,1,NULL,&da);
100:   DMSetOptionsPrefix(red,"da_");
101:   DMCompositeAddDM(packer,(DM)da);
102:   DMSetApplicationContext(packer,&user);

104:   packer->ops->creatematrix = DMCreateMatrix_MF;

106:   /* create nonlinear multi-level solver */
107:   SNESCreate(PETSC_COMM_WORLD,&snes);
108:   SNESSetDM(snes,packer);
109:   SNESSetFunction(snes,NULL,ComputeFunction,NULL);
110:   SNESSetJacobian(snes,NULL, NULL,ComputeJacobian_MF,NULL);

112:   SNESSetFromOptions(snes);

114:   if (use_monitor) {
115:     /* create graphics windows */
116:     PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u_lambda - state variables and Lagrange multipliers",-1,-1,-1,-1,&user.u_lambda_viewer);
117:     PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu_lambda - derivate w.r.t. state variables and Lagrange multipliers",-1,-1,-1,-1,&user.fu_lambda_viewer);
118:     SNESMonitorSet(snes,Monitor,0,0);
119:   }

121:   SNESSolve(snes,NULL,NULL);
122:   SNESDestroy(&snes);

124:   DMDestroy(&red);
125:   DMDestroy(&da);
126:   DMDestroy(&packer);
127:   if (use_monitor) {
128:     PetscViewerDestroy(&user.u_lambda_viewer);
129:     PetscViewerDestroy(&user.fu_lambda_viewer);
130:   }
131:   PetscFinalize();
132:   return 0;
133: }

135: typedef struct {
136:   PetscScalar u;
137:   PetscScalar lambda;
138: } ULambda;

142: /*
143:       Evaluates FU = Gradiant(L(w,u,lambda))

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

148: */
149: PetscErrorCode ComputeFunction(SNES snes,Vec U,Vec FU,void *ctx)
150: {
152:   PetscInt       xs,xm,i,N;
153:   ULambda        *u_lambda,*fu_lambda;
154:   PetscScalar    d,h,*w,*fw;
155:   Vec            vw,vfw,vu_lambda,vfu_lambda;
156:   DM             packer,red,da;

159:   VecGetDM(U, &packer);
160:   DMCompositeGetEntries(packer,&red,&da);
161:   DMCompositeGetLocalVectors(packer,&vw,&vu_lambda);
162:   DMCompositeScatter(packer,U,vw,vu_lambda);
163:   DMCompositeGetAccess(packer,FU,&vfw,&vfu_lambda);

165:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
166:   DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0);
167:   VecGetArray(vw,&w);
168:   VecGetArray(vfw,&fw);
169:   DMDAVecGetArray(da,vu_lambda,&u_lambda);
170:   DMDAVecGetArray(da,vfu_lambda,&fu_lambda);
171:   d    = N-1.0;
172:   h    = 1.0/d;

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

179:   /* derivative of L() w.r.t. u */
180:   for (i=xs; i<xs+xm; i++) {
181:     if      (i == 0)   fu_lambda[0].lambda   =    h*u_lambda[0].u   + 2.*d*u_lambda[0].lambda   - d*u_lambda[1].lambda;
182:     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;
183:     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;
184:     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;
185:     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);
186:   }

188:   /* derivative of L() w.r.t. lambda */
189:   for (i=xs; i<xs+xm; i++) {
190:     if      (i == 0)   fu_lambda[0].u   = 2.0*d*(u_lambda[0].u - w[0]);
191:     else if (i == N-1) fu_lambda[N-1].u = 2.0*d*u_lambda[N-1].u;
192:     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);
193:   }

195:   VecRestoreArray(vw,&w);
196:   VecRestoreArray(vfw,&fw);
197:   DMDAVecRestoreArray(da,vu_lambda,&u_lambda);
198:   DMDAVecRestoreArray(da,vfu_lambda,&fu_lambda);
199:   DMCompositeRestoreLocalVectors(packer,&vw,&vu_lambda);
200:   DMCompositeRestoreAccess(packer,FU,&vfw,&vfu_lambda);
201:   PetscLogFlops(13.0*N);
202:   return(0);
203: }

207: /*
208:     Computes the exact solution
209: */
210: PetscErrorCode u_solution(void *dummy,PetscInt n,const PetscScalar *x,PetscScalar *u)
211: {
212:   PetscInt i;

215:   for (i=0; i<n; i++) u[2*i] = x[i]*x[i] - 1.25*x[i] + .25;
216:   return(0);
217: }

221: PetscErrorCode ExactSolution(DM packer,Vec U)
222: {
223:   PF             pf;
224:   Vec            x,u_global;
225:   PetscScalar    *w;
226:   DM             da;
228:   PetscInt       m;

231:   DMCompositeGetEntries(packer,&m,&da);

233:   PFCreate(PETSC_COMM_WORLD,1,2,&pf);
234:   /* The cast through PETSC_UINTPTR_T is so that compilers will warn about casting to void * from void(*)(void) */
235:   PFSetType(pf,PFQUICK,(void*)(PETSC_UINTPTR_T)u_solution);
236:   DMGetCoordinates(da,&x);
237:   if (!x) {
238:     DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);
239:     DMGetCoordinates(da,&x);
240:   }
241:   DMCompositeGetAccess(packer,U,&w,&u_global,0);
242:   if (w) w[0] = .25;
243:   PFApplyVec(pf,x,u_global);
244:   PFDestroy(&pf);
245:   DMCompositeRestoreAccess(packer,U,&w,&u_global,0);
246:   return(0);
247: }

251: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy)
252: {
253:   UserCtx        *user;
255:   PetscInt       m,N;
256:   PetscScalar    *w,*dw;
257:   Vec            u_lambda,U,F,Uexact;
258:   DM             packer;
259:   PetscReal      norm;
260:   DM             da;

263:   SNESGetDM(snes,&packer);
264:   DMGetApplicationContext(packer,&user);
265:   SNESGetSolution(snes,&U);
266:   DMCompositeGetAccess(packer,U,&w,&u_lambda);
267:   VecView(u_lambda,user->u_lambda_viewer);
268:   DMCompositeRestoreAccess(packer,U,&w,&u_lambda);

270:   SNESGetFunction(snes,&F,0,0);
271:   DMCompositeGetAccess(packer,F,&w,&u_lambda);
272:   /* VecView(u_lambda,user->fu_lambda_viewer); */
273:   DMCompositeRestoreAccess(packer,U,&w,&u_lambda);

275:   DMCompositeGetEntries(packer,&m,&da);
276:   DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0);
277:   VecDuplicate(U,&Uexact);
278:   ExactSolution(packer,Uexact);
279:   VecAXPY(Uexact,-1.0,U);
280:   DMCompositeGetAccess(packer,Uexact,&dw,&u_lambda);
281:   VecStrideNorm(u_lambda,0,NORM_2,&norm);
282:   norm = norm/PetscSqrtReal((PetscReal)N-1.);
283:   if (dw) PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Error at x = 0 %g\n",(double)norm,(double)PetscRealPart(dw[0]));
284:   VecView(u_lambda,user->fu_lambda_viewer);
285:   DMCompositeRestoreAccess(packer,Uexact,&dw,&u_lambda);
286:   VecDestroy(&Uexact);
287:   return(0);
288: }

292: PetscErrorCode DMCreateMatrix_MF(DM packer,Mat *A)
293: {
295:   Vec            t;
296:   PetscInt       m;

299:   DMGetGlobalVector(packer,&t);
300:   VecGetLocalSize(t,&m);
301:   DMRestoreGlobalVector(packer,&t);
302:   MatCreateMFFD(PETSC_COMM_WORLD,m,m,PETSC_DETERMINE,PETSC_DETERMINE,A);
303:   MatSetUp(*A);
304:   return(0);
305: }

309: PetscErrorCode ComputeJacobian_MF(SNES snes,Vec x,Mat A,Mat B,void *ctx)
310: {

314:   MatMFFDSetFunction(A,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction,snes);
315:   MatMFFDSetBase(A,x,NULL);
316:   return(0);
317: }