Actual source code: ex21.c
petsc-3.4.5 2014-06-29
2: static const char help[] = "Solves PDE optimization problem using full-space method, treats state and adjoint variables separately.\n\n";
4: #include <petscdmda.h>
5: #include <petscdmredundant.h>
6: #include <petscdmcomposite.h>
7: #include <petscpf.h>
8: #include <petscsnes.h>
10: /*
12: w - design variables (what we change to get an optimal solution)
13: u - state variables (i.e. the PDE solution)
14: lambda - the Lagrange multipliers
16: U = (w u lambda)
18: fu, fw, flambda contain the gradient of L(w,u,lambda)
20: FU = (fw fu flambda)
22: In this example the PDE is
23: Uxx = 2,
24: u(0) = w(0), thus this is the free parameter
25: u(1) = 0
26: the function we wish to minimize is
27: \integral u^{2}
29: The exact solution for u is given by u(x) = x*x - 1.25*x + .25
31: Use the usual centered finite differences.
33: Note we treat the problem as non-linear though it happens to be linear
35: See ex22.c for the same code, but that interlaces the u and the lambda
37: */
39: typedef struct {
40: DM red1,da1,da2;
41: DM packer;
42: PetscViewer u_viewer,lambda_viewer;
43: PetscViewer fu_viewer,flambda_viewer;
44: } UserCtx;
46: extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
47: extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
52: int main(int argc,char **argv)
53: {
55: PetscInt its;
56: Vec U,FU;
57: SNES snes;
58: UserCtx user;
60: PetscInitialize(&argc,&argv,(char*)0,help);
62: /* Create a global vector that includes a single redundant array and two da arrays */
63: DMCompositeCreate(PETSC_COMM_WORLD,&user.packer);
64: DMRedundantCreate(PETSC_COMM_WORLD,0,1,&user.red1);
65: DMCompositeAddDM(user.packer,user.red1);
66: DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-5,1,1,NULL,&user.da1);
67: DMCompositeAddDM(user.packer,user.da1);
68: DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-5,1,1,NULL,&user.da2);
69: DMCompositeAddDM(user.packer,user.da2);
70: DMCreateGlobalVector(user.packer,&U);
71: VecDuplicate(U,&FU);
73: /* create graphics windows */
74: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u - state variables",-1,-1,-1,-1,&user.u_viewer);
75: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"lambda - Lagrange multipliers",-1,-1,-1,-1,&user.lambda_viewer);
76: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu - derivate w.r.t. state variables",-1,-1,-1,-1,&user.fu_viewer);
77: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"flambda - derivate w.r.t. Lagrange multipliers",-1,-1,-1,-1,&user.flambda_viewer);
80: /* create nonlinear solver */
81: SNESCreate(PETSC_COMM_WORLD,&snes);
82: SNESSetFunction(snes,FU,FormFunction,&user);
83: SNESSetFromOptions(snes);
84: SNESMonitorSet(snes,Monitor,&user,0);
85: SNESSolve(snes,NULL,U);
86: SNESGetIterationNumber(snes,&its);
87: SNESDestroy(&snes);
89: DMDestroy(&user.red1);
90: DMDestroy(&user.da1);
91: DMDestroy(&user.da2);
92: DMDestroy(&user.packer);
93: VecDestroy(&U);
94: VecDestroy(&FU);
95: PetscViewerDestroy(&user.u_viewer);
96: PetscViewerDestroy(&user.lambda_viewer);
97: PetscViewerDestroy(&user.fu_viewer);
98: PetscViewerDestroy(&user.flambda_viewer);
99: PetscFinalize();
100: return 0;
101: }
105: /*
106: Evaluates FU = Gradiant(L(w,u,lambda))
108: */
109: PetscErrorCode FormFunction(SNES snes,Vec U,Vec FU,void *dummy)
110: {
111: UserCtx *user = (UserCtx*)dummy;
113: PetscInt xs,xm,i,N;
114: PetscScalar *u,*lambda,*w,*fu,*fw,*flambda,d,h;
115: Vec vw,vu,vlambda,vfw,vfu,vflambda;
118: DMCompositeGetLocalVectors(user->packer,&vw,&vu,&vlambda);
119: DMCompositeGetLocalVectors(user->packer,&vfw,&vfu,&vflambda);
120: DMCompositeScatter(user->packer,U,vw,vu,vlambda);
122: DMDAGetCorners(user->da1,&xs,NULL,NULL,&xm,NULL,NULL);
123: DMDAGetInfo(user->da1,0,&N,0,0,0,0,0,0,0,0,0,0,0);
124: VecGetArray(vw,&w);
125: VecGetArray(vfw,&fw);
126: DMDAVecGetArray(user->da1,vu,&u);
127: DMDAVecGetArray(user->da1,vfu,&fu);
128: DMDAVecGetArray(user->da1,vlambda,&lambda);
129: DMDAVecGetArray(user->da1,vflambda,&flambda);
130: d = (N-1.0);
131: h = 1.0/d;
133: /* derivative of L() w.r.t. w */
134: if (xs == 0) { /* only first processor computes this */
135: fw[0] = -2.*d*lambda[0];
136: }
138: /* derivative of L() w.r.t. u */
139: for (i=xs; i<xs+xm; i++) {
140: if (i == 0) flambda[0] = h*u[0] + 2.*d*lambda[0] - d*lambda[1];
141: else if (i == 1) flambda[1] = 2.*h*u[1] + 2.*d*lambda[1] - d*lambda[2];
142: else if (i == N-1) flambda[N-1] = h*u[N-1] + 2.*d*lambda[N-1] - d*lambda[N-2];
143: else if (i == N-2) flambda[N-2] = 2.*h*u[N-2] + 2.*d*lambda[N-2] - d*lambda[N-3];
144: else flambda[i] = 2.*h*u[i] - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]);
145: }
147: /* derivative of L() w.r.t. lambda */
148: for (i=xs; i<xs+xm; i++) {
149: if (i == 0) fu[0] = 2.0*d*(u[0] - w[0]);
150: else if (i == N-1) fu[N-1] = 2.0*d*u[N-1];
151: else fu[i] = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h);
152: }
154: VecRestoreArray(vw,&w);
155: VecRestoreArray(vfw,&fw);
156: DMDAVecRestoreArray(user->da1,vu,&u);
157: DMDAVecRestoreArray(user->da1,vfu,&fu);
158: DMDAVecRestoreArray(user->da1,vlambda,&lambda);
159: DMDAVecRestoreArray(user->da1,vflambda,&flambda);
161: DMCompositeGather(user->packer,FU,INSERT_VALUES,vfw,vfu,vflambda);
162: DMCompositeRestoreLocalVectors(user->packer,&vw,&vu,&vlambda);
163: DMCompositeRestoreLocalVectors(user->packer,&vfw,&vfu,&vflambda);
164: return(0);
165: }
169: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy)
170: {
171: UserCtx *user = (UserCtx*)dummy;
173: Vec w,u,lambda,U,F;
176: SNESGetSolution(snes,&U);
177: DMCompositeGetAccess(user->packer,U,&w,&u,&lambda);
178: VecView(u,user->u_viewer);
179: VecView(lambda,user->lambda_viewer);
180: DMCompositeRestoreAccess(user->packer,U,&w,&u,&lambda);
182: SNESGetFunction(snes,&F,0,0);
183: DMCompositeGetAccess(user->packer,F,&w,&u,&lambda);
184: VecView(u,user->fu_viewer);
185: VecView(lambda,user->flambda_viewer);
186: DMCompositeRestoreAccess(user->packer,F,&w,&u,&lambda);
187: return(0);
188: }