Actual source code: ex21.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2: static const char help[] = "Solves PDE optimization problem using full-space method, treats state and adjoint variables separately.\n\n";

  4:  #include <petscdm.h>
  5:  #include <petscdmda.h>
  6:  #include <petscdmredundant.h>
  7:  #include <petscdmcomposite.h>
  8:  #include <petscpf.h>
  9:  #include <petscsnes.h>

 11: /*

 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

 17:             U = (w u lambda)

 19:        fu, fw, flambda contain the gradient of L(w,u,lambda)

 21:             FU = (fw fu flambda)

 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}

 30:        The exact solution for u is given by u(x) = x*x - 1.25*x + .25

 32:        Use the usual centered finite differences.

 34:        Note we treat the problem as non-linear though it happens to be linear

 36:        See ex22.c for the same code, but that interlaces the u and the lambda

 38: */

 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;

 47: extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
 48: extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);


 51: int main(int argc,char **argv)
 52: {
 54:   PetscInt       its;
 55:   Vec            U,FU;
 56:   SNES           snes;
 57:   UserCtx        user;

 59:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 61:   /* Create a global vector that includes a single redundant array and two da arrays */
 62:   DMCompositeCreate(PETSC_COMM_WORLD,&user.packer);
 63:   DMRedundantCreate(PETSC_COMM_WORLD,0,1,&user.red1);
 64:   DMCompositeAddDM(user.packer,user.red1);
 65:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,1,1,NULL,&user.da1);
 66:   DMSetFromOptions(user.da1);
 67:   DMSetUp(user.da1);
 68:   DMCompositeAddDM(user.packer,user.da1);
 69:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,1,1,NULL,&user.da2);
 70:   DMSetFromOptions(user.da2);
 71:   DMSetUp(user.da2);
 72:   DMDASetFieldName(user.da1,0,"u");
 73:   DMDASetFieldName(user.da2,0,"lambda");
 74:   DMCompositeAddDM(user.packer,user.da2);
 75:   DMCreateGlobalVector(user.packer,&U);
 76:   VecDuplicate(U,&FU);

 78:   /* create graphics windows */
 79:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u - state variables",-1,-1,-1,-1,&user.u_viewer);
 80:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"lambda - Lagrange multipliers",-1,-1,-1,-1,&user.lambda_viewer);
 81:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu - derivate w.r.t. state variables",-1,-1,-1,-1,&user.fu_viewer);
 82:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"flambda - derivate w.r.t. Lagrange multipliers",-1,-1,-1,-1,&user.flambda_viewer);


 85:   /* create nonlinear solver */
 86:   SNESCreate(PETSC_COMM_WORLD,&snes);
 87:   SNESSetFunction(snes,FU,FormFunction,&user);
 88:   SNESSetFromOptions(snes);
 89:   SNESMonitorSet(snes,Monitor,&user,0);
 90:   SNESSolve(snes,NULL,U);
 91:   SNESGetIterationNumber(snes,&its);
 92:   SNESDestroy(&snes);

 94:   DMDestroy(&user.red1);
 95:   DMDestroy(&user.da1);
 96:   DMDestroy(&user.da2);
 97:   DMDestroy(&user.packer);
 98:   VecDestroy(&U);
 99:   VecDestroy(&FU);
100:   PetscViewerDestroy(&user.u_viewer);
101:   PetscViewerDestroy(&user.lambda_viewer);
102:   PetscViewerDestroy(&user.fu_viewer);
103:   PetscViewerDestroy(&user.flambda_viewer);
104:   PetscFinalize();
105:   return ierr;
106: }

108: /*
109:       Evaluates FU = Gradiant(L(w,u,lambda))

111: */
112: PetscErrorCode FormFunction(SNES snes,Vec U,Vec FU,void *dummy)
113: {
114:   UserCtx        *user = (UserCtx*)dummy;
116:   PetscInt       xs,xm,i,N;
117:   PetscScalar    *u,*lambda,*w,*fu,*fw,*flambda,d,h;
118:   Vec            vw,vu,vlambda,vfw,vfu,vflambda;

121:   DMCompositeGetLocalVectors(user->packer,&vw,&vu,&vlambda);
122:   DMCompositeGetLocalVectors(user->packer,&vfw,&vfu,&vflambda);
123:   DMCompositeScatter(user->packer,U,vw,vu,vlambda);

125:   DMDAGetCorners(user->da1,&xs,NULL,NULL,&xm,NULL,NULL);
126:   DMDAGetInfo(user->da1,0,&N,0,0,0,0,0,0,0,0,0,0,0);
127:   VecGetArray(vw,&w);
128:   VecGetArray(vfw,&fw);
129:   DMDAVecGetArray(user->da1,vu,&u);
130:   DMDAVecGetArray(user->da1,vfu,&fu);
131:   DMDAVecGetArray(user->da1,vlambda,&lambda);
132:   DMDAVecGetArray(user->da1,vflambda,&flambda);
133:   d    = (N-1.0);
134:   h    = 1.0/d;

136:   /* derivative of L() w.r.t. w */
137:   if (xs == 0) { /* only first processor computes this */
138:     fw[0] = -2.*d*lambda[0];
139:   }

141:   /* derivative of L() w.r.t. u */
142:   for (i=xs; i<xs+xm; i++) {
143:     if      (i == 0)   flambda[0]   =    h*u[0]   + 2.*d*lambda[0]   - d*lambda[1];
144:     else if (i == 1)   flambda[1]   = 2.*h*u[1]   + 2.*d*lambda[1]   - d*lambda[2];
145:     else if (i == N-1) flambda[N-1] =    h*u[N-1] + 2.*d*lambda[N-1] - d*lambda[N-2];
146:     else if (i == N-2) flambda[N-2] = 2.*h*u[N-2] + 2.*d*lambda[N-2] - d*lambda[N-3];
147:     else               flambda[i]   = 2.*h*u[i]   - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]);
148:   }

150:   /* derivative of L() w.r.t. lambda */
151:   for (i=xs; i<xs+xm; i++) {
152:     if      (i == 0)   fu[0]   = 2.0*d*(u[0] - w[0]);
153:     else if (i == N-1) fu[N-1] = 2.0*d*u[N-1];
154:     else               fu[i]   = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h);
155:   }

157:   VecRestoreArray(vw,&w);
158:   VecRestoreArray(vfw,&fw);
159:   DMDAVecRestoreArray(user->da1,vu,&u);
160:   DMDAVecRestoreArray(user->da1,vfu,&fu);
161:   DMDAVecRestoreArray(user->da1,vlambda,&lambda);
162:   DMDAVecRestoreArray(user->da1,vflambda,&flambda);

164:   DMCompositeGather(user->packer,INSERT_VALUES,FU,vfw,vfu,vflambda);
165:   DMCompositeRestoreLocalVectors(user->packer,&vw,&vu,&vlambda);
166:   DMCompositeRestoreLocalVectors(user->packer,&vfw,&vfu,&vflambda);
167:   return(0);
168: }

170: PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy)
171: {
172:   UserCtx        *user = (UserCtx*)dummy;
174:   Vec            w,u,lambda,U,F;

177:   SNESGetSolution(snes,&U);
178:   DMCompositeGetAccess(user->packer,U,&w,&u,&lambda);
179:   VecView(u,user->u_viewer);
180:   VecView(lambda,user->lambda_viewer);
181:   DMCompositeRestoreAccess(user->packer,U,&w,&u,&lambda);

183:   SNESGetFunction(snes,&F,0,0);
184:   DMCompositeGetAccess(user->packer,F,&w,&u,&lambda);
185:   VecView(u,user->fu_viewer);
186:   VecView(lambda,user->flambda_viewer);
187:   DMCompositeRestoreAccess(user->packer,F,&w,&u,&lambda);
188:   return(0);
189: }

191: /*TEST

193:    test:
194:       nsize: 4
195:       args: -snes_linesearch_monitor -snes_mf -snes_monitor_short -nox -ksp_monitor_short -snes_converged_reason
196:       requires: !single

198: TEST*/