xref: /petsc/src/snes/tutorials/ex21.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static const char help[] = "Solves PDE optimization problem using full-space method, treats state and adjoint variables separately.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown #include <petscdmredundant.h>
7c4762a1bSJed Brown #include <petscdmcomposite.h>
8c4762a1bSJed Brown #include <petscpf.h>
9c4762a1bSJed Brown #include <petscsnes.h>
10c4762a1bSJed Brown 
11c4762a1bSJed Brown /*
12c4762a1bSJed Brown 
13c4762a1bSJed Brown        w - design variables (what we change to get an optimal solution)
14c4762a1bSJed Brown        u - state variables (i.e. the PDE solution)
15c4762a1bSJed Brown        lambda - the Lagrange multipliers
16c4762a1bSJed Brown 
17c4762a1bSJed Brown             U = (w u lambda)
18c4762a1bSJed Brown 
19c4762a1bSJed Brown        fu, fw, flambda contain the gradient of L(w,u,lambda)
20c4762a1bSJed Brown 
21c4762a1bSJed Brown             FU = (fw fu flambda)
22c4762a1bSJed Brown 
23c4762a1bSJed Brown        In this example the PDE is
24c4762a1bSJed Brown                              Uxx = 2,
25c4762a1bSJed Brown                             u(0) = w(0), thus this is the free parameter
26c4762a1bSJed Brown                             u(1) = 0
27c4762a1bSJed Brown        the function we wish to minimize is
28c4762a1bSJed Brown                             \integral u^{2}
29c4762a1bSJed Brown 
30c4762a1bSJed Brown        The exact solution for u is given by u(x) = x*x - 1.25*x + .25
31c4762a1bSJed Brown 
32c4762a1bSJed Brown        Use the usual centered finite differences.
33c4762a1bSJed Brown 
34c4762a1bSJed Brown        Note we treat the problem as non-linear though it happens to be linear
35c4762a1bSJed Brown 
36c4762a1bSJed Brown        See ex22.c for the same code, but that interlaces the u and the lambda
37c4762a1bSJed Brown 
38c4762a1bSJed Brown */
39c4762a1bSJed Brown 
40c4762a1bSJed Brown typedef struct {
41c4762a1bSJed Brown   DM          red1,da1,da2;
42c4762a1bSJed Brown   DM          packer;
43c4762a1bSJed Brown   PetscViewer u_viewer,lambda_viewer;
44c4762a1bSJed Brown   PetscViewer fu_viewer,flambda_viewer;
45c4762a1bSJed Brown } UserCtx;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
48c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
49c4762a1bSJed Brown 
50c4762a1bSJed Brown int main(int argc,char **argv)
51c4762a1bSJed Brown {
52c4762a1bSJed Brown   PetscInt       its;
53c4762a1bSJed Brown   Vec            U,FU;
54c4762a1bSJed Brown   SNES           snes;
55c4762a1bSJed Brown   UserCtx        user;
56c4762a1bSJed Brown 
57*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Create a global vector that includes a single redundant array and two da arrays */
605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeCreate(PETSC_COMM_WORLD,&user.packer));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRedundantCreate(PETSC_COMM_WORLD,0,1,&user.red1));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeAddDM(user.packer,user.red1));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,1,1,NULL,&user.da1));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(user.da1));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(user.da1));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeAddDM(user.packer,user.da1));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,1,1,NULL,&user.da2));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(user.da2));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(user.da2));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(user.da1,0,"u"));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(user.da2,0,"lambda"));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeAddDM(user.packer,user.da2));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(user.packer,&U));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(U,&FU));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* create graphics windows */
775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u - state variables",-1,-1,-1,-1,&user.u_viewer));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"lambda - Lagrange multipliers",-1,-1,-1,-1,&user.lambda_viewer));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu - derivate w.r.t. state variables",-1,-1,-1,-1,&user.fu_viewer));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"flambda - derivate w.r.t. Lagrange multipliers",-1,-1,-1,-1,&user.flambda_viewer));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* create nonlinear solver */
835f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,FU,FormFunction,&user));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESMonitorSet(snes,Monitor,&user,0));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,U));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
90c4762a1bSJed Brown 
915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&user.red1));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&user.da1));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&user.da2));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&user.packer));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&FU));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&user.u_viewer));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&user.lambda_viewer));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&user.fu_viewer));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&user.flambda_viewer));
101*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
102*b122ec5aSJacob Faibussowitsch   return 0;
103c4762a1bSJed Brown }
104c4762a1bSJed Brown 
105c4762a1bSJed Brown /*
106c4762a1bSJed Brown       Evaluates FU = Gradiant(L(w,u,lambda))
107c4762a1bSJed Brown 
108c4762a1bSJed Brown */
109c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec U,Vec FU,void *dummy)
110c4762a1bSJed Brown {
111c4762a1bSJed Brown   UserCtx        *user = (UserCtx*)dummy;
112c4762a1bSJed Brown   PetscInt       xs,xm,i,N;
113c4762a1bSJed Brown   PetscScalar    *u,*lambda,*w,*fu,*fw,*flambda,d,h;
114c4762a1bSJed Brown   Vec            vw,vu,vlambda,vfw,vfu,vflambda;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   PetscFunctionBeginUser;
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeGetLocalVectors(user->packer,&vw,&vu,&vlambda));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeGetLocalVectors(user->packer,&vfw,&vfu,&vflambda));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeScatter(user->packer,U,vw,vu,vlambda));
120c4762a1bSJed Brown 
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user->da1,&xs,NULL,NULL,&xm,NULL,NULL));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(user->da1,0,&N,0,0,0,0,0,0,0,0,0,0,0));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(vw,&w));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(vfw,&fw));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(user->da1,vu,&u));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(user->da1,vfu,&fu));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(user->da1,vlambda,&lambda));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(user->da1,vflambda,&flambda));
129c4762a1bSJed Brown   d    = (N-1.0);
130c4762a1bSJed Brown   h    = 1.0/d;
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* derivative of L() w.r.t. w */
133c4762a1bSJed Brown   if (xs == 0) { /* only first processor computes this */
134c4762a1bSJed Brown     fw[0] = -2.*d*lambda[0];
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* derivative of L() w.r.t. u */
138c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
139c4762a1bSJed Brown     if      (i == 0)   flambda[0]   =    h*u[0]   + 2.*d*lambda[0]   - d*lambda[1];
140c4762a1bSJed Brown     else if (i == 1)   flambda[1]   = 2.*h*u[1]   + 2.*d*lambda[1]   - d*lambda[2];
141c4762a1bSJed Brown     else if (i == N-1) flambda[N-1] =    h*u[N-1] + 2.*d*lambda[N-1] - d*lambda[N-2];
142c4762a1bSJed Brown     else if (i == N-2) flambda[N-2] = 2.*h*u[N-2] + 2.*d*lambda[N-2] - d*lambda[N-3];
143c4762a1bSJed Brown     else               flambda[i]   = 2.*h*u[i]   - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]);
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* derivative of L() w.r.t. lambda */
147c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
148c4762a1bSJed Brown     if      (i == 0)   fu[0]   = 2.0*d*(u[0] - w[0]);
149c4762a1bSJed Brown     else if (i == N-1) fu[N-1] = 2.0*d*u[N-1];
150c4762a1bSJed Brown     else               fu[i]   = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h);
151c4762a1bSJed Brown   }
152c4762a1bSJed Brown 
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(vw,&w));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(vfw,&fw));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(user->da1,vu,&u));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(user->da1,vfu,&fu));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(user->da1,vlambda,&lambda));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(user->da1,vflambda,&flambda));
159c4762a1bSJed Brown 
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeGather(user->packer,INSERT_VALUES,FU,vfw,vfu,vflambda));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeRestoreLocalVectors(user->packer,&vw,&vu,&vlambda));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeRestoreLocalVectors(user->packer,&vfw,&vfu,&vflambda));
163c4762a1bSJed Brown   PetscFunctionReturn(0);
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy)
167c4762a1bSJed Brown {
168c4762a1bSJed Brown   UserCtx        *user = (UserCtx*)dummy;
169c4762a1bSJed Brown   Vec            w,u,lambda,U,F;
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   PetscFunctionBeginUser;
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetSolution(snes,&U));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeGetAccess(user->packer,U,&w,&u,&lambda));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(u,user->u_viewer));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(lambda,user->lambda_viewer));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeRestoreAccess(user->packer,U,&w,&u,&lambda));
177c4762a1bSJed Brown 
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetFunction(snes,&F,0,0));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeGetAccess(user->packer,F,&w,&u,&lambda));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(u,user->fu_viewer));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(lambda,user->flambda_viewer));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCompositeRestoreAccess(user->packer,F,&w,&u,&lambda));
183c4762a1bSJed Brown   PetscFunctionReturn(0);
184c4762a1bSJed Brown }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown /*TEST
187c4762a1bSJed Brown 
188c4762a1bSJed Brown    test:
189c4762a1bSJed Brown       nsize: 4
190c4762a1bSJed Brown       args: -snes_linesearch_monitor -snes_mf -pc_type none -snes_monitor_short -nox -ksp_monitor_short -snes_converged_reason
191c4762a1bSJed Brown       requires: !single
192c4762a1bSJed Brown 
193c4762a1bSJed Brown TEST*/
194