Actual source code: ex21.c
1: /*$Id: ex21.c,v 1.11 2001/08/07 03:04:16 balay Exp $*/
3: static char help[] = "Solves PDE optimization problem.nn";
5: #include petscda.h
6: #include petscpf.h
7: #include petscsnes.h
9: /*
11: w - design variables (what we change to get an optimal solution)
12: u - state variables (i.e. the PDE solution)
13: lambda - the Lagrange multipliers
15: U = (w u lambda)
17: fu, fw, flambda contain the gradient of L(w,u,lambda)
19: FU = (fw fu flambda)
21: In this example the PDE is
22: Uxx = 2,
23: u(0) = w(0), thus this is the free parameter
24: u(1) = 0
25: the function we wish to minimize is
26: integral u^{2}
28: The exact solution for u is given by u(x) = x*x - 1.25*x + .25
30: Use the usual centered finite differences.
32: Note we treat the problem as non-linear though it happens to be linear
34: See ex22.c for the same code, but that interlaces the u and the lambda
36: */
38: typedef struct {
39: DA da1,da2;
40: int nredundant;
41: VecPack packer;
42: PetscViewer u_viewer,lambda_viewer;
43: PetscViewer fu_viewer,flambda_viewer;
44: } UserCtx;
46: extern int FormFunction(SNES,Vec,Vec,void*);
47: extern int Monitor(SNES,int,PetscReal,void*);
50: int main(int argc,char **argv)
51: {
52: int ierr,its;
53: Vec U,FU;
54: SNES snes;
55: UserCtx user;
57: PetscInitialize(&argc,&argv,(char*)0,help);
59: /* Create a global vector that includes a single redundant array and two da arrays */
60: VecPackCreate(PETSC_COMM_WORLD,&user.packer);
61: user.nredundant = 1;
62: VecPackAddArray(user.packer,user.nredundant);
63: DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,-5,1,1,PETSC_NULL,&user.da1);
64: VecPackAddDA(user.packer,user.da1);
65: DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,-5,1,1,PETSC_NULL,&user.da2);
66: VecPackAddDA(user.packer,user.da2);
67: VecPackCreateGlobalVector(user.packer,&U);
68: VecDuplicate(U,&FU);
70: /* create graphics windows */
71: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u - state variables",-1,-1,-1,-1,&user.u_viewer);
72: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"lambda - Lagrange multipliers",-1,-1,-1,-1,&user.lambda_viewer);
73: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu - derivate w.r.t. state variables",-1,-1,-1,-1,&user.fu_viewer);
74: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"flambda - derivate w.r.t. Lagrange multipliers",-1,-1,-1,-1,&user.flambda_viewer);
77: /* create nonlinear solver */
78: SNESCreate(PETSC_COMM_WORLD,&snes);
79: SNESSetFunction(snes,FU,FormFunction,&user);
80: SNESSetFromOptions(snes);
81: SNESSetMonitor(snes,Monitor,&user,0);
82: SNESSolve(snes,U,&its);
83: SNESDestroy(snes);
85: DADestroy(user.da1);
86: DADestroy(user.da2);
87: VecPackDestroy(user.packer);
88: VecDestroy(U);
89: VecDestroy(FU);
90: PetscViewerDestroy(user.u_viewer);
91: PetscViewerDestroy(user.lambda_viewer);
92: PetscViewerDestroy(user.fu_viewer);
93: PetscViewerDestroy(user.flambda_viewer);
94: PetscFinalize();
95: return 0;
96: }
97:
98: /*
99: Evaluates FU = Gradiant(L(w,u,lambda))
101: */
102: int FormFunction(SNES snes,Vec U,Vec FU,void* dummy)
103: {
104: UserCtx *user = (UserCtx*)dummy;
105: int ierr,xs,xm,i,N;
106: PetscScalar *u,*lambda,*w,*fu,*fw,*flambda,d,h;
107: Vec vu,vlambda,vfu,vflambda;
110: VecPackGetLocalVectors(user->packer,&w,&vu,&vlambda);
111: VecPackGetLocalVectors(user->packer,&fw,&vfu,&vflambda);
112: VecPackScatter(user->packer,U,w,vu,vlambda);
114: DAGetCorners(user->da1,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
115: DAGetInfo(user->da1,0,&N,0,0,0,0,0,0,0,0,0);
116: DAVecGetArray(user->da1,vu,(void**)&u);
117: DAVecGetArray(user->da1,vfu,(void**)&fu);
118: DAVecGetArray(user->da1,vlambda,(void**)&lambda);
119: DAVecGetArray(user->da1,vflambda,(void**)&flambda);
120: d = (N-1.0);
121: h = 1.0/d;
123: /* derivative of L() w.r.t. w */
124: if (xs == 0) { /* only first processor computes this */
125: fw[0] = -2.*d*lambda[0];
126: }
128: /* derivative of L() w.r.t. u */
129: for (i=xs; i<xs+xm; i++) {
130: if (i == 0) flambda[0] = h*u[0] + 2.*d*lambda[0] - d*lambda[1];
131: else if (i == 1) flambda[1] = 2.*h*u[1] + 2.*d*lambda[1] - d*lambda[2];
132: else if (i == N-1) flambda[N-1] = h*u[N-1] + 2.*d*lambda[N-1] - d*lambda[N-2];
133: else if (i == N-2) flambda[N-2] = 2.*h*u[N-2] + 2.*d*lambda[N-2] - d*lambda[N-3];
134: else flambda[i] = 2.*h*u[i] - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]);
135: }
137: /* derivative of L() w.r.t. lambda */
138: for (i=xs; i<xs+xm; i++) {
139: if (i == 0) fu[0] = 2.0*d*(u[0] - w[0]);
140: else if (i == N-1) fu[N-1] = 2.0*d*u[N-1];
141: else fu[i] = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h);
142: }
144: DAVecRestoreArray(user->da1,vu,(void**)&u);
145: DAVecRestoreArray(user->da1,vfu,(void**)&fu);
146: DAVecRestoreArray(user->da1,vlambda,(void**)&lambda);
147: DAVecRestoreArray(user->da1,vflambda,(void**)&flambda);
149: VecPackGather(user->packer,FU,fw,vfu,vflambda);
150: VecPackRestoreLocalVectors(user->packer,&w,&vu,&vlambda);
151: VecPackRestoreLocalVectors(user->packer,&fw,&vfu,&vflambda);
152: return(0);
153: }
155: int Monitor(SNES snes,int its,PetscReal rnorm,void *dummy)
156: {
157: UserCtx *user = (UserCtx*)dummy;
158: int ierr;
159: PetscScalar *w;
160: Vec u,lambda,U,F;
163: SNESGetSolution(snes,&U);
164: VecPackGetAccess(user->packer,U,&w,&u,&lambda);
165: VecView(u,user->u_viewer);
166: VecView(lambda,user->lambda_viewer);
167: VecPackRestoreAccess(user->packer,U,&w,&u,&lambda);
169: SNESGetFunction(snes,&F,0,0);
170: VecPackGetAccess(user->packer,F,&w,&u,&lambda);
171: VecView(u,user->fu_viewer);
172: VecView(lambda,user->flambda_viewer);
173: VecPackRestoreAccess(user->packer,F,&w,&u,&lambda);
174: return(0);
175: }