Actual source code: ex24.c
1: /*$Id: ex24.c,v 1.15 2001/04/10 19:37:05 bsmith Exp $*/
3: static char help[] = "Solves PDE optimization problem of ex22.c with finite differences for adjoint.nn";
5: #include petscda.h
6: #include petscpf.h
7: #include petscsnes.h
9: /*
11: Minimize F(w,u) such that G(w,u) = 0
13: L(w,u,lambda) = F(w,u) + lambda^T G(w,u)
15: w - design variables (what we change to get an optimal solution)
16: u - state variables (i.e. the PDE solution)
17: lambda - the Lagrange multipliers
19: U = (w u lambda)
21: fu, fw, flambda contain the gradient of L(w,u,lambda)
23: FU = (fw fu flambda)
25: In this example the PDE is
26: Uxx = 2,
27: u(0) = w(0), thus this is the free parameter
28: u(1) = 0
29: the function we wish to minimize is
30: integral u^{2}
32: The exact solution for u is given by u(x) = x*x - 1.25*x + .25
34: Use the usual centered finite differences.
36: Note we treat the problem as non-linear though it happens to be linear
38: The lambda and u are NOT interlaced.
39: */
41: extern int FormFunction(SNES,Vec,Vec,void*);
42: extern int PDEFormFunction(Scalar*,Vec,Vec,DA);
44: int main(int argc,char **argv)
45: {
46: int ierr,N = 5,nlevels,i;
47: DA da;
48: DMMG *dmmg;
49: VecPack packer;
50: Mat J;
53: PetscInitialize(&argc,&argv,PETSC_NULL,help);
54: PetscOptionsGetInt(PETSC_NULL,"-N",&N,PETSC_NULL);
56: /* Hardwire several options; can be changed at command line */
57: PetscOptionsSetValue("-dmmg_grid_sequence",PETSC_NULL);
58: PetscOptionsSetValue("-ksp_type","fgmres");
59: PetscOptionsSetValue("-ksp_max_it","5");
60: PetscOptionsSetValue("-pc_mg_type","full");
61: PetscOptionsSetValue("-mg_coarse_ksp_type","gmres");
62: PetscOptionsSetValue("-mg_levels_ksp_type","gmres");
63: PetscOptionsSetValue("-mg_coarse_ksp_max_it","6");
64: PetscOptionsSetValue("-mg_levels_ksp_max_it","3");
65: PetscOptionsSetValue("-snes_mf_type","wp");
66: PetscOptionsSetValue("-snes_mf_compute_norma","no");
67: PetscOptionsSetValue("-snes_mf_compute_normu","no");
68: PetscOptionsSetValue("-snes_eq_ls","basic");
69: /* PetscOptionsSetValue("-snes_eq_ls","basicnonorms"); */
70: PetscOptionsInsert(&argc,&argv,PETSC_NULL);
71:
72: /* Create a global vector from a da arrays */
73: DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,N,1,1,PETSC_NULL,&da);
74: VecPackCreate(PETSC_COMM_WORLD,&packer);
75: VecPackAddArray(packer,1);
76: VecPackAddDA(packer,da);
77: VecPackAddDA(packer,da);
78: DADestroy(da);
80: /* create nonlinear multi-level solver */
81: DMMGCreate(PETSC_COMM_WORLD,2,PETSC_NULL,&dmmg);
82: DMMGSetUseMatrixFree(dmmg);
83: DMMGSetDM(dmmg,(DM)packer);
84: VecPackDestroy(packer);
86: /* Create Jacobian of PDE function for each level */
87: nlevels = DMMGGetLevels(dmmg);
88: for (i=0; i<nlevels; i++) {
89: packer = (VecPack)dmmg[i]->dm;
90: ierr = VecPackGetEntries(packer,PETSC_NULL,&da,PETSC_NULL);
91: ierr = DAGetColoring(da,IS_COLORING_GLOBAL,MATMPIAIJ,PETSC_NULL,&J);
92: dmmg[i]->user = (void*)J;
93: }
95: DMMGSetSNES(dmmg,FormFunction,PETSC_NULL);
96: DMMGSolve(dmmg);
98: VecView(DMMGGetx(dmmg),PETSC_VIEWER_SOCKET_WORLD);
99: for (i=0; i<nlevels; i++) {
100: MatDestroy((Mat)dmmg[i]->user);
101: }
102: DMMGDestroy(dmmg);
104: PetscFinalize();
105: return 0;
106: }
107:
108: /*
109: Enforces the PDE on the grid
110: This local function acts on the ghosted version of U (accessed via DAGetLocalVector())
111: BUT the global, nonghosted version of FU
113: */
114: int PDEFormFunction(Scalar *w,Vec vu,Vec vfu,DA da)
115: {
116: int ierr,xs,xm,i,N;
117: Scalar *u,*fu,d,h;
121: DAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
122: DAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0);
123: DAVecGetArray(da,vu,(void**)&u);
124: DAVecGetArray(da,vfu,(void**)&fu);
125: d = N-1.0;
126: h = 1.0/d;
128: for (i=xs; i<xs+xm; i++) {
129: if (i == 0) fu[0] = 2.0*d*(u[0] - w[0]) + h*u[0]*u[0];
130: else if (i == N-1) fu[N-1] = 2.0*d*u[N-1] + h*u[N-1]*u[N-1];
131: else fu[i] = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h) + h*u[i]*u[i];
132: }
134: DAVecRestoreArray(da,vu,(void**)&u);
135: DAVecRestoreArray(da,vfu,(void**)&fu);
136: PetscLogFlops(9*N);
137: return(0);
138: }
141: /*
142: Evaluates FU = Gradiant(L(w,u,lambda))
144: This local function acts on the ghosted version of U (accessed via VecPackGetLocalVectors() and
145: VecPackScatter()) BUT the global, nonghosted version of FU (via VecPackAccess()).
147: This function uses PDEFormFunction() to enforce the PDE constraint equations and its adjoint
148: for the Lagrange multiplier equations
150: */
151: int FormFunction(SNES snes,Vec U,Vec FU,void* dummy)
152: {
153: DMMG dmmg = (DMMG)dummy;
154: int ierr,xs,xm,i,N,nredundant;
155: Scalar *u,*w,*fw,*fu,*lambda,*flambda,d,h,h2;
156: Vec vu,vlambda,vfu,vflambda,vglambda;
157: DA da,dadummy;
158: VecPack packer = (VecPack)dmmg->dm;
161: VecPackGetEntries(packer,&nredundant,&da,&dadummy);
162: VecPackGetLocalVectors(packer,&w,&vu,&vlambda);
163: VecPackScatter(packer,U,w,vu,vlambda);
164: VecPackGetAccess(packer,FU,&fw,&vfu,&vflambda);
165: VecPackGetAccess(packer,U,0,0,&vglambda);
167: /* Evaluate the Jacobian of PDEFormFunction() */
168: /* PDEFormJacobian(da,w,vu,(Mat)dmmg->user);
169: MatMultTranspose((Mat)dmmg->user,vglambda,vflambda); */
171: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_MATLAB);
172: /* MatView((Mat)dmmg->user,PETSC_VIEWER_STDOUT_WORLD); */
173: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
175: /* derivative of constraint portion of L() w.r.t. u */
176: PDEFormFunction(w,vu,vfu,da);
178: DAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
179: DAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0);
180: DAVecGetArray(da,vu,(void**)&u);
181: DAVecGetArray(da,vfu,(void**)&fu);
182: DAVecGetArray(da,vlambda,(void**)&lambda);
183: DAVecGetArray(da,vflambda,(void**)&flambda);
184: d = (N-1.0);
185: h = 1.0/d;
186: h2 = 2.0*h;
188: /* derivative of L() w.r.t. w */
189: if (xs == 0) { /* only first processor computes this */
190: fw[0] = -2.*d*lambda[0];
191: }
193: for (i=xs; i<xs+xm; i++) {
194: if (i == 0) flambda[0] = 2.*d*lambda[0] - d*lambda[1] + h2*lambda[0]*u[0];
195: else if (i == 1) flambda[1] = 2.*d*lambda[1] - d*lambda[2] + h2*lambda[1]*u[1];
196: else if (i == N-1) flambda[N-1] = 2.*d*lambda[N-1] - d*lambda[N-2] + h2*lambda[N-1]*u[N-1];
197: else if (i == N-2) flambda[N-2] = 2.*d*lambda[N-2] - d*lambda[N-3] + h2*lambda[N-2]*u[N-2];
198: else flambda[i] = - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]) + h2*lambda[i]*u[i];
199: }
201: /* derivative of function part of L() w.r.t. u */
202: for (i=xs; i<xs+xm; i++) {
203: if (i == 0) flambda[0] += h*u[0];
204: else if (i == 1) flambda[1] += h2*u[1];
205: else if (i == N-1) flambda[N-1] += h*u[N-1];
206: else if (i == N-2) flambda[N-2] += h2*u[N-2];
207: else flambda[i] += h2*u[i];
208: }
210: DAVecRestoreArray(da,vu,(void**)&u);
211: DAVecRestoreArray(da,vfu,(void**)&fu);
212: DAVecRestoreArray(da,vlambda,(void**)&lambda);
213: DAVecRestoreArray(da,vflambda,(void**)&flambda);
216: VecPackRestoreLocalVectors(packer,&w,&vu,&vlambda);
217: VecPackRestoreAccess(packer,FU,&fw,&vfu,&vflambda);
218: VecPackRestoreAccess(packer,U,0,0,&vglambda);
220: PetscLogFlops(9*N);
221: return(0);
222: }