xref: /petsc/src/snes/tutorials/ex21.c (revision 81b34b0469c4c918ac4feb79d232c68b446ecfbd)
1 
2 static const char help[] = "Solves PDE optimization problem using full-space method, treats state and adjoint variables separately.\n\n";
3 
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 
11 /*
12 
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
16 
17             U = (w u lambda)
18 
19        fu, fw, flambda contain the gradient of L(w,u,lambda)
20 
21             FU = (fw fu flambda)
22 
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}
29 
30        The exact solution for u is given by u(x) = x*x - 1.25*x + .25
31 
32        Use the usual centered finite differences.
33 
34        Note we treat the problem as non-linear though it happens to be linear
35 
36        See ex22.c for the same code, but that interlaces the u and the lambda
37 
38 */
39 
40 typedef struct {
41   DM          red1,da1,da2;
42   DM          packer;
43   PetscViewer u_viewer,lambda_viewer;
44   PetscViewer fu_viewer,flambda_viewer;
45 } UserCtx;
46 
47 extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
48 extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
49 
50 int main(int argc,char **argv)
51 {
52   PetscErrorCode ierr;
53   PetscInt       its;
54   Vec            U,FU;
55   SNES           snes;
56   UserCtx        user;
57 
58   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
59 
60   /* Create a global vector that includes a single redundant array and two da arrays */
61   ierr = DMCompositeCreate(PETSC_COMM_WORLD,&user.packer);CHKERRQ(ierr);
62   ierr = DMRedundantCreate(PETSC_COMM_WORLD,0,1,&user.red1);CHKERRQ(ierr);
63   ierr = DMCompositeAddDM(user.packer,user.red1);CHKERRQ(ierr);
64   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,1,1,NULL,&user.da1);CHKERRQ(ierr);
65   ierr = DMSetFromOptions(user.da1);CHKERRQ(ierr);
66   ierr = DMSetUp(user.da1);CHKERRQ(ierr);
67   ierr = DMCompositeAddDM(user.packer,user.da1);CHKERRQ(ierr);
68   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,1,1,NULL,&user.da2);CHKERRQ(ierr);
69   ierr = DMSetFromOptions(user.da2);CHKERRQ(ierr);
70   ierr = DMSetUp(user.da2);CHKERRQ(ierr);
71   ierr = DMDASetFieldName(user.da1,0,"u");CHKERRQ(ierr);
72   ierr = DMDASetFieldName(user.da2,0,"lambda");CHKERRQ(ierr);
73   ierr = DMCompositeAddDM(user.packer,user.da2);CHKERRQ(ierr);
74   ierr = DMCreateGlobalVector(user.packer,&U);CHKERRQ(ierr);
75   ierr = VecDuplicate(U,&FU);CHKERRQ(ierr);
76 
77   /* create graphics windows */
78   ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u - state variables",-1,-1,-1,-1,&user.u_viewer);CHKERRQ(ierr);
79   ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"lambda - Lagrange multipliers",-1,-1,-1,-1,&user.lambda_viewer);CHKERRQ(ierr);
80   ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu - derivate w.r.t. state variables",-1,-1,-1,-1,&user.fu_viewer);CHKERRQ(ierr);
81   ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"flambda - derivate w.r.t. Lagrange multipliers",-1,-1,-1,-1,&user.flambda_viewer);CHKERRQ(ierr);
82 
83   /* create nonlinear solver */
84   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
85   ierr = SNESSetFunction(snes,FU,FormFunction,&user);CHKERRQ(ierr);
86   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
87   ierr = SNESMonitorSet(snes,Monitor,&user,0);CHKERRQ(ierr);
88   ierr = SNESSolve(snes,NULL,U);CHKERRQ(ierr);
89   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
90   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
91 
92   ierr = DMDestroy(&user.red1);CHKERRQ(ierr);
93   ierr = DMDestroy(&user.da1);CHKERRQ(ierr);
94   ierr = DMDestroy(&user.da2);CHKERRQ(ierr);
95   ierr = DMDestroy(&user.packer);CHKERRQ(ierr);
96   ierr = VecDestroy(&U);CHKERRQ(ierr);
97   ierr = VecDestroy(&FU);CHKERRQ(ierr);
98   ierr = PetscViewerDestroy(&user.u_viewer);CHKERRQ(ierr);
99   ierr = PetscViewerDestroy(&user.lambda_viewer);CHKERRQ(ierr);
100   ierr = PetscViewerDestroy(&user.fu_viewer);CHKERRQ(ierr);
101   ierr = PetscViewerDestroy(&user.flambda_viewer);CHKERRQ(ierr);
102   ierr = PetscFinalize();
103   return ierr;
104 }
105 
106 /*
107       Evaluates FU = Gradiant(L(w,u,lambda))
108 
109 */
110 PetscErrorCode FormFunction(SNES snes,Vec U,Vec FU,void *dummy)
111 {
112   UserCtx        *user = (UserCtx*)dummy;
113   PetscErrorCode ierr;
114   PetscInt       xs,xm,i,N;
115   PetscScalar    *u,*lambda,*w,*fu,*fw,*flambda,d,h;
116   Vec            vw,vu,vlambda,vfw,vfu,vflambda;
117 
118   PetscFunctionBeginUser;
119   ierr = DMCompositeGetLocalVectors(user->packer,&vw,&vu,&vlambda);CHKERRQ(ierr);
120   ierr = DMCompositeGetLocalVectors(user->packer,&vfw,&vfu,&vflambda);CHKERRQ(ierr);
121   ierr = DMCompositeScatter(user->packer,U,vw,vu,vlambda);CHKERRQ(ierr);
122 
123   ierr = DMDAGetCorners(user->da1,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
124   ierr = DMDAGetInfo(user->da1,0,&N,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
125   ierr = VecGetArray(vw,&w);CHKERRQ(ierr);
126   ierr = VecGetArray(vfw,&fw);CHKERRQ(ierr);
127   ierr = DMDAVecGetArray(user->da1,vu,&u);CHKERRQ(ierr);
128   ierr = DMDAVecGetArray(user->da1,vfu,&fu);CHKERRQ(ierr);
129   ierr = DMDAVecGetArray(user->da1,vlambda,&lambda);CHKERRQ(ierr);
130   ierr = DMDAVecGetArray(user->da1,vflambda,&flambda);CHKERRQ(ierr);
131   d    = (N-1.0);
132   h    = 1.0/d;
133 
134   /* derivative of L() w.r.t. w */
135   if (xs == 0) { /* only first processor computes this */
136     fw[0] = -2.*d*lambda[0];
137   }
138 
139   /* derivative of L() w.r.t. u */
140   for (i=xs; i<xs+xm; i++) {
141     if      (i == 0)   flambda[0]   =    h*u[0]   + 2.*d*lambda[0]   - d*lambda[1];
142     else if (i == 1)   flambda[1]   = 2.*h*u[1]   + 2.*d*lambda[1]   - d*lambda[2];
143     else if (i == N-1) flambda[N-1] =    h*u[N-1] + 2.*d*lambda[N-1] - d*lambda[N-2];
144     else if (i == N-2) flambda[N-2] = 2.*h*u[N-2] + 2.*d*lambda[N-2] - d*lambda[N-3];
145     else               flambda[i]   = 2.*h*u[i]   - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]);
146   }
147 
148   /* derivative of L() w.r.t. lambda */
149   for (i=xs; i<xs+xm; i++) {
150     if      (i == 0)   fu[0]   = 2.0*d*(u[0] - w[0]);
151     else if (i == N-1) fu[N-1] = 2.0*d*u[N-1];
152     else               fu[i]   = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h);
153   }
154 
155   ierr = VecRestoreArray(vw,&w);CHKERRQ(ierr);
156   ierr = VecRestoreArray(vfw,&fw);CHKERRQ(ierr);
157   ierr = DMDAVecRestoreArray(user->da1,vu,&u);CHKERRQ(ierr);
158   ierr = DMDAVecRestoreArray(user->da1,vfu,&fu);CHKERRQ(ierr);
159   ierr = DMDAVecRestoreArray(user->da1,vlambda,&lambda);CHKERRQ(ierr);
160   ierr = DMDAVecRestoreArray(user->da1,vflambda,&flambda);CHKERRQ(ierr);
161 
162   ierr = DMCompositeGather(user->packer,INSERT_VALUES,FU,vfw,vfu,vflambda);CHKERRQ(ierr);
163   ierr = DMCompositeRestoreLocalVectors(user->packer,&vw,&vu,&vlambda);CHKERRQ(ierr);
164   ierr = DMCompositeRestoreLocalVectors(user->packer,&vfw,&vfu,&vflambda);CHKERRQ(ierr);
165   PetscFunctionReturn(0);
166 }
167 
168 PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy)
169 {
170   UserCtx        *user = (UserCtx*)dummy;
171   PetscErrorCode ierr;
172   Vec            w,u,lambda,U,F;
173 
174   PetscFunctionBeginUser;
175   ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr);
176   ierr = DMCompositeGetAccess(user->packer,U,&w,&u,&lambda);CHKERRQ(ierr);
177   ierr = VecView(u,user->u_viewer);CHKERRQ(ierr);
178   ierr = VecView(lambda,user->lambda_viewer);CHKERRQ(ierr);
179   ierr = DMCompositeRestoreAccess(user->packer,U,&w,&u,&lambda);CHKERRQ(ierr);
180 
181   ierr = SNESGetFunction(snes,&F,0,0);CHKERRQ(ierr);
182   ierr = DMCompositeGetAccess(user->packer,F,&w,&u,&lambda);CHKERRQ(ierr);
183   ierr = VecView(u,user->fu_viewer);CHKERRQ(ierr);
184   ierr = VecView(lambda,user->flambda_viewer);CHKERRQ(ierr);
185   ierr = DMCompositeRestoreAccess(user->packer,F,&w,&u,&lambda);CHKERRQ(ierr);
186   PetscFunctionReturn(0);
187 }
188 
189 /*TEST
190 
191    test:
192       nsize: 4
193       args: -snes_linesearch_monitor -snes_mf -pc_type none -snes_monitor_short -nox -ksp_monitor_short -snes_converged_reason
194       requires: !single
195 
196 TEST*/
197