Actual source code: ex22.c
petsc-3.3-p7 2013-05-11
2: static const char help[] = "Solves PDE optimization problem using full-space method, interlaces state and adjoint variables.\n\n";
4: #include <petscdmda.h>
5: #include <petscdmredundant.h>
6: #include <petscdmcomposite.h>
7: #include <petscpf.h>
8: #include <petscsnes.h>
9: #include <petsc-private/dmimpl.h>
11: /*
13: w - design variables (what we change to get an optimal solution)
14: u - state variables (i.e. the PDE solution)
15: lambda - the Lagrange multipliers
17: U = (w [u_0 lambda_0 u_1 lambda_1 .....])
19: fu, fw, flambda contain the gradient of L(w,u,lambda)
21: FU = (fw [fu_0 flambda_0 .....])
23: In this example the PDE is
24: Uxx = 2,
25: u(0) = w(0), thus this is the free parameter
26: u(1) = 0
27: the function we wish to minimize is
28: \integral u^{2}
30: The exact solution for u is given by u(x) = x*x - 1.25*x + .25
32: Use the usual centered finite differences.
34: Note we treat the problem as non-linear though it happens to be linear
36: See ex21.c for the same code, but that does NOT interlaces the u and the lambda
38: The vectors u_lambda and fu_lambda contain the u and the lambda interlaced
39: */
41: typedef struct {
42: PetscViewer u_lambda_viewer;
43: PetscViewer fu_lambda_viewer;
44: } UserCtx;
46: extern PetscErrorCode ComputeFunction(SNES,Vec,Vec,void*);
47: extern PetscErrorCode ComputeJacobian_MF(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
48: extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
50: /*
51: Uses full multigrid preconditioner with GMRES (with no preconditioner inside the GMRES) as the
52: smoother on all levels. This is because (1) in the matrix free case no matrix entries are
53: available for doing Jacobi or SOR preconditioning and (2) the explicit matrix case the diagonal
54: entry for the control variable is zero which means default SOR will not work.
56: */
57: char common_options[] = "-ksp_type fgmres\
58: -snes_grid_sequence 4 \
59: -pc_type mg\
60: -mg_levels_pc_type none \
61: -mg_coarse_pc_type none \
62: -pc_mg_type full \
63: -mg_coarse_ksp_type gmres \
64: -mg_levels_ksp_type gmres \
65: -mg_coarse_ksp_max_it 6 \
66: -mg_levels_ksp_max_it 3";
68: char matrix_free_options[] = "-mat_mffd_compute_normu no \
69: -mat_mffd_type wp";
71: extern PetscErrorCode DMCreateMatrix_MF(DM,const MatType,Mat*);
75: int main(int argc,char **argv)
76: {
78: UserCtx user;
79: DM red,da;
80: SNES snes;
81: DM packer;
82: PetscBool use_monitor = PETSC_FALSE;
84: PetscInitialize(&argc,&argv,PETSC_NULL,help);
85: PetscOptionsSetFromOptions();
87: /* Hardwire several options; can be changed at command line */
88: PetscOptionsInsertString(common_options);
89: PetscOptionsInsertString(matrix_free_options);
90: PetscOptionsInsert(&argc,&argv,PETSC_NULL);
91: PetscOptionsGetBool(PETSC_NULL,"-use_monitor",&use_monitor,PETSC_IGNORE);
93: /* Create a global vector that includes a single redundant array and two da arrays */
94: DMCompositeCreate(PETSC_COMM_WORLD,&packer);
95: DMRedundantCreate(PETSC_COMM_WORLD,0,1,&red);
96: DMSetOptionsPrefix(red,"red_");
97: DMCompositeAddDM(packer,red);
98: DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-5,2,1,PETSC_NULL,&da);
99: DMSetOptionsPrefix(red,"da_");
100: DMCompositeAddDM(packer,(DM)da);
101: DMSetApplicationContext(packer,&user);
102: DMSNESSetFunction(packer,ComputeFunction,PETSC_NULL);
103: DMSNESSetJacobian(packer,ComputeJacobian_MF,PETSC_NULL);
104: packer->ops->creatematrix = DMCreateMatrix_MF;
106: /* create nonlinear multi-level solver */
107: SNESCreate(PETSC_COMM_WORLD,&snes);
108: SNESSetDM(snes,packer);
109: SNESSetFromOptions(snes);
111: if (use_monitor) {
112: /* create graphics windows */
113: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u_lambda - state variables and Lagrange multipliers",-1,-1,-1,-1,&user.u_lambda_viewer);
114: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu_lambda - derivate w.r.t. state variables and Lagrange multipliers",-1,-1,-1,-1,&user.fu_lambda_viewer);
115: SNESMonitorSet(snes,Monitor,0,0);
116: }
118: SNESSolve(snes,PETSC_NULL,PETSC_NULL);
119: SNESDestroy(&snes);
121: DMDestroy(&red);
122: DMDestroy(&da);
123: DMDestroy(&packer);
124: if (use_monitor) {
125: PetscViewerDestroy(&user.u_lambda_viewer);
126: PetscViewerDestroy(&user.fu_lambda_viewer);
127: }
128: PetscFinalize();
129: return 0;
130: }
132: typedef struct {
133: PetscScalar u;
134: PetscScalar lambda;
135: } ULambda;
139: /*
140: Evaluates FU = Gradiant(L(w,u,lambda))
142: This local function acts on the ghosted version of U (accessed via DMCompositeGetLocalVectors() and
143: DMCompositeScatter()) BUT the global, nonghosted version of FU (via DMCompositeGetAccess()).
145: */
146: PetscErrorCode ComputeFunction(SNES snes,Vec U,Vec FU,void *ctx)
147: {
149: PetscInt xs,xm,i,N;
150: ULambda *u_lambda,*fu_lambda;
151: PetscScalar d,h,*w,*fw;
152: Vec vw,vfw,vu_lambda,vfu_lambda;
153: DM packer,red,da;
156: PetscObjectQuery((PetscObject)U,"DM",(PetscObject*)&packer); /* Ugly way to get context */
157: DMCompositeGetEntries(packer,&red,&da);
158: DMCompositeGetLocalVectors(packer,&vw,&vu_lambda);
159: DMCompositeScatter(packer,U,vw,vu_lambda);
160: DMCompositeGetAccess(packer,FU,&vfw,&vfu_lambda);
162: DMDAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
163: DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0);
164: VecGetArray(vw,&w);
165: VecGetArray(vfw,&fw);
166: DMDAVecGetArray(da,vu_lambda,&u_lambda);
167: DMDAVecGetArray(da,vfu_lambda,&fu_lambda);
168: d = N-1.0;
169: h = 1.0/d;
171: /* derivative of L() w.r.t. w */
172: if (xs == 0) { /* only first processor computes this */
173: fw[0] = -2.0*d*u_lambda[0].lambda;
174: }
176: /* derivative of L() w.r.t. u */
177: for (i=xs; i<xs+xm; i++) {
178: if (i == 0) fu_lambda[0].lambda = h*u_lambda[0].u + 2.*d*u_lambda[0].lambda - d*u_lambda[1].lambda;
179: 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;
180: 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;
181: 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;
182: 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);
183: }
185: /* derivative of L() w.r.t. lambda */
186: for (i=xs; i<xs+xm; i++) {
187: if (i == 0) fu_lambda[0].u = 2.0*d*(u_lambda[0].u - w[0]);
188: else if (i == N-1) fu_lambda[N-1].u = 2.0*d*u_lambda[N-1].u;
189: 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);
190: }
192: VecRestoreArray(vw,&w);
193: VecRestoreArray(vfw,&fw);
194: DMDAVecRestoreArray(da,vu_lambda,&u_lambda);
195: DMDAVecRestoreArray(da,vfu_lambda,&fu_lambda);
196: DMCompositeRestoreLocalVectors(packer,&vw,&vu_lambda);
197: DMCompositeRestoreAccess(packer,FU,&vfw,&vfu_lambda);
198: PetscLogFlops(13.0*N);
199: return(0);
200: }
204: /*
205: Computes the exact solution
206: */
207: PetscErrorCode u_solution(void *dummy,PetscInt n,const PetscScalar *x,PetscScalar *u)
208: {
209: PetscInt i;
211: for (i=0; i<n; i++) {
212: u[2*i] = x[i]*x[i] - 1.25*x[i] + .25;
213: }
214: return(0);
215: }
219: PetscErrorCode ExactSolution(DM packer,Vec U)
220: {
221: PF pf;
222: Vec x,u_global;
223: PetscScalar *w;
224: DM da;
226: PetscInt m;
229: DMCompositeGetEntries(packer,&m,&da);
231: PFCreate(PETSC_COMM_WORLD,1,2,&pf);
232: PFSetType(pf,PFQUICK,(void*)u_solution);
233: DMDAGetCoordinates(da,&x);
234: if (!x) {
235: DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);
236: DMDAGetCoordinates(da,&x);
237: }
238: DMCompositeGetAccess(packer,U,&w,&u_global,0);
239: if (w) w[0] = .25;
240: PFApplyVec(pf,x,u_global);
241: PFDestroy(&pf);
242: DMCompositeRestoreAccess(packer,U,&w,&u_global,0);
243: return(0);
244: }
248: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy)
249: {
250: UserCtx *user;
252: PetscInt m,N;
253: PetscScalar *w,*dw;
254: Vec u_lambda,U,F,Uexact;
255: DM packer;
256: PetscReal norm;
257: DM da;
260: SNESGetDM(snes,&packer);
261: DMGetApplicationContext(packer,&user);
262: SNESGetSolution(snes,&U);
263: DMCompositeGetAccess(packer,U,&w,&u_lambda);
264: VecView(u_lambda,user->u_lambda_viewer);
265: DMCompositeRestoreAccess(packer,U,&w,&u_lambda);
267: SNESGetFunction(snes,&F,0,0);
268: DMCompositeGetAccess(packer,F,&w,&u_lambda);
269: /* VecView(u_lambda,user->fu_lambda_viewer); */
270: DMCompositeRestoreAccess(packer,U,&w,&u_lambda);
272: DMCompositeGetEntries(packer,&m,&da);
273: DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0);
274: VecDuplicate(U,&Uexact);
275: ExactSolution(packer,Uexact);
276: VecAXPY(Uexact,-1.0,U);
277: DMCompositeGetAccess(packer,Uexact,&dw,&u_lambda);
278: VecStrideNorm(u_lambda,0,NORM_2,&norm);
279: norm = norm/sqrt(N-1.);
280: if (dw) PetscPrintf(PETSC_COMM_WORLD,"Norm of error %G Error at x = 0 %G\n",norm,PetscRealPart(dw[0]));
281: VecView(u_lambda,user->fu_lambda_viewer);
282: DMCompositeRestoreAccess(packer,Uexact,&dw,&u_lambda);
283: VecDestroy(&Uexact);
284: return(0);
285: }
289: PetscErrorCode DMCreateMatrix_MF(DM packer,const MatType stype,Mat *A)
290: {
292: Vec t;
293: PetscInt m;
296: DMGetGlobalVector(packer,&t);
297: VecGetLocalSize(t,&m);
298: DMRestoreGlobalVector(packer,&t);
299: MatCreateMFFD(PETSC_COMM_WORLD,m,m,PETSC_DETERMINE,PETSC_DETERMINE,A);
300: MatSetUp(*A);
301: return(0);
302: }
306: PetscErrorCode ComputeJacobian_MF(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,void *ctx)
307: {
311: MatMFFDSetFunction(*A,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction,snes);
312: MatMFFDSetBase(*A,x,PETSC_NULL);
313: return(0);
314: }