Actual source code: ex22.c
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*);
74: int main(int argc,char **argv)
75: {
77: UserCtx user;
78: DM red,da;
79: SNES snes;
80: DM packer;
81: PetscBool use_monitor = PETSC_FALSE;
83: PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
85: /* Hardwire several options; can be changed at command line */
86: PetscOptionsInsertString(NULL,common_options);
87: PetscOptionsInsertString(NULL,matrix_free_options);
88: PetscOptionsInsert(NULL,&argc,&argv,NULL);
89: PetscOptionsGetBool(NULL,NULL,"-use_monitor",&use_monitor,PETSC_IGNORE);
91: /* Create a global vector that includes a single redundant array and two da arrays */
92: DMCompositeCreate(PETSC_COMM_WORLD,&packer);
93: DMRedundantCreate(PETSC_COMM_WORLD,0,1,&red);
94: DMSetOptionsPrefix(red,"red_");
95: DMCompositeAddDM(packer,red);
96: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,2,1,NULL,&da);
97: DMSetOptionsPrefix(red,"da_");
98: DMSetFromOptions(da);
99: DMSetUp(da);
100: DMDASetFieldName(da,0,"u");
101: DMDASetFieldName(da,1,"lambda");
102: DMCompositeAddDM(packer,(DM)da);
103: DMSetApplicationContext(packer,&user);
105: packer->ops->creatematrix = DMCreateMatrix_MF;
107: /* create nonlinear multi-level solver */
108: SNESCreate(PETSC_COMM_WORLD,&snes);
109: SNESSetDM(snes,packer);
110: SNESSetFunction(snes,NULL,ComputeFunction,NULL);
111: SNESSetJacobian(snes,NULL, NULL,ComputeJacobian_MF,NULL);
113: SNESSetFromOptions(snes);
115: if (use_monitor) {
116: /* create graphics windows */
117: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u_lambda - state variables and Lagrange multipliers",-1,-1,-1,-1,&user.u_lambda_viewer);
118: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu_lambda - derivate w.r.t. state variables and Lagrange multipliers",-1,-1,-1,-1,&user.fu_lambda_viewer);
119: SNESMonitorSet(snes,Monitor,0,0);
120: }
122: SNESSolve(snes,NULL,NULL);
123: SNESDestroy(&snes);
125: DMDestroy(&red);
126: DMDestroy(&da);
127: DMDestroy(&packer);
128: if (use_monitor) {
129: PetscViewerDestroy(&user.u_lambda_viewer);
130: PetscViewerDestroy(&user.fu_lambda_viewer);
131: }
132: PetscFinalize();
133: return ierr;
134: }
136: typedef struct {
137: PetscScalar u;
138: PetscScalar lambda;
139: } ULambda;
141: /*
142: Evaluates FU = Gradiant(L(w,u,lambda))
144: This local function acts on the ghosted version of U (accessed via DMCompositeGetLocalVectors() and
145: DMCompositeScatter()) BUT the global, nonghosted version of FU (via DMCompositeGetAccess()).
147: */
148: PetscErrorCode ComputeFunction(SNES snes,Vec U,Vec FU,void *ctx)
149: {
151: PetscInt xs,xm,i,N;
152: ULambda *u_lambda,*fu_lambda;
153: PetscScalar d,h,*w,*fw;
154: Vec vw,vfw,vu_lambda,vfu_lambda;
155: DM packer,red,da;
158: VecGetDM(U, &packer);
159: DMCompositeGetEntries(packer,&red,&da);
160: DMCompositeGetLocalVectors(packer,&vw,&vu_lambda);
161: DMCompositeScatter(packer,U,vw,vu_lambda);
162: DMCompositeGetAccess(packer,FU,&vfw,&vfu_lambda);
164: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
165: DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0);
166: VecGetArray(vw,&w);
167: VecGetArray(vfw,&fw);
168: DMDAVecGetArray(da,vu_lambda,&u_lambda);
169: DMDAVecGetArray(da,vfu_lambda,&fu_lambda);
170: d = N-1.0;
171: h = 1.0/d;
173: /* derivative of L() w.r.t. w */
174: if (xs == 0) { /* only first processor computes this */
175: fw[0] = -2.0*d*u_lambda[0].lambda;
176: }
178: /* derivative of L() w.r.t. u */
179: for (i=xs; i<xs+xm; i++) {
180: if (i == 0) fu_lambda[0].lambda = h*u_lambda[0].u + 2.*d*u_lambda[0].lambda - d*u_lambda[1].lambda;
181: 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;
182: 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;
183: 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;
184: 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);
185: }
187: /* derivative of L() w.r.t. lambda */
188: for (i=xs; i<xs+xm; i++) {
189: if (i == 0) fu_lambda[0].u = 2.0*d*(u_lambda[0].u - w[0]);
190: else if (i == N-1) fu_lambda[N-1].u = 2.0*d*u_lambda[N-1].u;
191: 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);
192: }
194: VecRestoreArray(vw,&w);
195: VecRestoreArray(vfw,&fw);
196: DMDAVecRestoreArray(da,vu_lambda,&u_lambda);
197: DMDAVecRestoreArray(da,vfu_lambda,&fu_lambda);
198: DMCompositeRestoreLocalVectors(packer,&vw,&vu_lambda);
199: DMCompositeRestoreAccess(packer,FU,&vfw,&vfu_lambda);
200: PetscLogFlops(13.0*N);
201: return(0);
202: }
204: /*
205: Computes the exact solution
206: */
207: PetscErrorCode u_solution(void *dummy,PetscInt n,const PetscScalar *x,PetscScalar *u)
208: {
209: PetscInt i;
212: for (i=0; i<n; i++) u[2*i] = x[i]*x[i] - 1.25*x[i] + .25;
213: return(0);
214: }
216: PetscErrorCode ExactSolution(DM packer,Vec U)
217: {
218: PF pf;
219: Vec x,u_global;
220: PetscScalar *w;
221: DM da;
223: PetscInt m;
226: DMCompositeGetEntries(packer,&m,&da);
228: PFCreate(PETSC_COMM_WORLD,1,2,&pf);
229: /* The cast through PETSC_UINTPTR_T is so that compilers will warn about casting to void * from void(*)(void) */
230: PFSetType(pf,PFQUICK,(void*)(PETSC_UINTPTR_T)u_solution);
231: DMGetCoordinates(da,&x);
232: if (!x) {
233: DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);
234: DMGetCoordinates(da,&x);
235: }
236: DMCompositeGetAccess(packer,U,&w,&u_global,0);
237: if (w) w[0] = .25;
238: PFApplyVec(pf,x,u_global);
239: PFDestroy(&pf);
240: DMCompositeRestoreAccess(packer,U,&w,&u_global,0);
241: return(0);
242: }
244: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy)
245: {
246: UserCtx *user;
248: PetscInt m,N;
249: PetscScalar *w,*dw;
250: Vec u_lambda,U,F,Uexact;
251: DM packer;
252: PetscReal norm;
253: DM da;
256: SNESGetDM(snes,&packer);
257: DMGetApplicationContext(packer,&user);
258: SNESGetSolution(snes,&U);
259: DMCompositeGetAccess(packer,U,&w,&u_lambda);
260: VecView(u_lambda,user->u_lambda_viewer);
261: DMCompositeRestoreAccess(packer,U,&w,&u_lambda);
263: SNESGetFunction(snes,&F,0,0);
264: DMCompositeGetAccess(packer,F,&w,&u_lambda);
265: /* VecView(u_lambda,user->fu_lambda_viewer); */
266: DMCompositeRestoreAccess(packer,U,&w,&u_lambda);
268: DMCompositeGetEntries(packer,&m,&da);
269: DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0);
270: VecDuplicate(U,&Uexact);
271: ExactSolution(packer,Uexact);
272: VecAXPY(Uexact,-1.0,U);
273: DMCompositeGetAccess(packer,Uexact,&dw,&u_lambda);
274: VecStrideNorm(u_lambda,0,NORM_2,&norm);
275: norm = norm/PetscSqrtReal((PetscReal)N-1.);
276: if (dw) {PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Error at x = 0 %g\n",(double)norm,(double)PetscRealPart(dw[0]));}
277: VecView(u_lambda,user->fu_lambda_viewer);
278: DMCompositeRestoreAccess(packer,Uexact,&dw,&u_lambda);
279: VecDestroy(&Uexact);
280: return(0);
281: }
283: PetscErrorCode DMCreateMatrix_MF(DM packer,Mat *A)
284: {
286: Vec t;
287: PetscInt m;
290: DMGetGlobalVector(packer,&t);
291: VecGetLocalSize(t,&m);
292: DMRestoreGlobalVector(packer,&t);
293: MatCreateMFFD(PETSC_COMM_WORLD,m,m,PETSC_DETERMINE,PETSC_DETERMINE,A);
294: MatSetUp(*A);
295: MatSetDM(*A,packer);
296: return(0);
297: }
299: PetscErrorCode ComputeJacobian_MF(SNES snes,Vec x,Mat A,Mat B,void *ctx)
300: {
304: MatMFFDSetFunction(A,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction,snes);
305: MatMFFDSetBase(A,x,NULL);
306: return(0);
307: }
309: /*TEST
311: test:
312: nsize: 2
313: args: -da_grid_x 10 -snes_converged_reason -ksp_converged_reason -snes_view
314: requires: !single
316: TEST*/