Actual source code: ex22.c

petsc-master 2016-12-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 2 \
 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);if (ierr) return ierr;
 86:   PetscOptionsSetFromOptions(NULL);

 88:   /* Hardwire several options; can be changed at command line */
 89:   PetscOptionsInsertString(NULL,common_options);
 90:   PetscOptionsInsertString(NULL,matrix_free_options);
 91:   PetscOptionsInsert(NULL,&argc,&argv,NULL);
 92:   PetscOptionsGetBool(NULL,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:   DMSetFromOptions(da);
102:   DMSetUp(da);
103:   DMCompositeAddDM(packer,(DM)da);
104:   DMSetApplicationContext(packer,&user);

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

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

114:   SNESSetFromOptions(snes);

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

123:   SNESSolve(snes,NULL,NULL);
124:   SNESDestroy(&snes);

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

137: typedef struct {
138:   PetscScalar u;
139:   PetscScalar lambda;
140: } ULambda;

144: /*
145:       Evaluates FU = Gradiant(L(w,u,lambda))

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

294: PetscErrorCode DMCreateMatrix_MF(DM packer,Mat *A)
295: {
297:   Vec            t;
298:   PetscInt       m;

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

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

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