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: #undef __FUNCT__
 52: int main(int argc,char **argv)
 53: {
 54:   int     ierr,its;
 55:   Vec     U,FU;
 56:   SNES    snes;
 57:   UserCtx user;

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

 61:   /* Create a global vector that includes a single redundant array and two da arrays */
 62:   VecPackCreate(PETSC_COMM_WORLD,&user.packer);
 63:   user.nredundant = 1;
 64:   VecPackAddArray(user.packer,user.nredundant);
 65:   DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,-5,1,1,PETSC_NULL,&user.da1);
 66:   VecPackAddDA(user.packer,user.da1);
 67:   DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,-5,1,1,PETSC_NULL,&user.da2);
 68:   VecPackAddDA(user.packer,user.da2);
 69:   VecPackCreateGlobalVector(user.packer,&U);
 70:   VecDuplicate(U,&FU);

 72:   /* create graphics windows */
 73:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u - state variables",-1,-1,-1,-1,&user.u_viewer);
 74:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"lambda - Lagrange multipliers",-1,-1,-1,-1,&user.lambda_viewer);
 75:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu - derivate w.r.t. state variables",-1,-1,-1,-1,&user.fu_viewer);
 76:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"flambda - derivate w.r.t. Lagrange multipliers",-1,-1,-1,-1,&user.flambda_viewer);


 79:   /* create nonlinear solver */
 80:   SNESCreate(PETSC_COMM_WORLD,&snes);
 81:   SNESSetFunction(snes,FU,FormFunction,&user);
 82:   SNESSetFromOptions(snes);
 83:   SNESSetMonitor(snes,Monitor,&user,0);
 84:   SNESSolve(snes,U,&its);
 85:   SNESDestroy(snes);

 87:   DADestroy(user.da1);
 88:   DADestroy(user.da2);
 89:   VecPackDestroy(user.packer);
 90:   VecDestroy(U);
 91:   VecDestroy(FU);
 92:   PetscViewerDestroy(user.u_viewer);
 93:   PetscViewerDestroy(user.lambda_viewer);
 94:   PetscViewerDestroy(user.fu_viewer);
 95:   PetscViewerDestroy(user.flambda_viewer);
 96:   PetscFinalize();
 97:   return 0;
 98: }
 99: 
100: /*
101:       Evaluates FU = Gradiant(L(w,u,lambda))

103: */
104: int FormFunction(SNES snes,Vec U,Vec FU,void* dummy)
105: {
106:   UserCtx *user = (UserCtx*)dummy;
107:   int     ierr,xs,xm,i,N;
108:   PetscScalar  *u,*lambda,*w,*fu,*fw,*flambda,d,h;
109:   Vec     vu,vlambda,vfu,vflambda;

112:   VecPackGetLocalVectors(user->packer,&w,&vu,&vlambda);
113:   VecPackGetLocalVectors(user->packer,&fw,&vfu,&vflambda);
114:   VecPackScatter(user->packer,U,w,vu,vlambda);

116:   DAGetCorners(user->da1,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
117:   DAGetInfo(user->da1,0,&N,0,0,0,0,0,0,0,0,0);
118:   DAVecGetArray(user->da1,vu,(void**)&u);
119:   DAVecGetArray(user->da1,vfu,(void**)&fu);
120:   DAVecGetArray(user->da1,vlambda,(void**)&lambda);
121:   DAVecGetArray(user->da1,vflambda,(void**)&flambda);
122:   d    = (N-1.0);
123:   h    = 1.0/d;

125:   /* derivative of L() w.r.t. w */
126:   if (xs == 0) { /* only first processor computes this */
127:     fw[0] = -2.*d*lambda[0];
128:   }

130:   /* derivative of L() w.r.t. u */
131:   for (i=xs; i<xs+xm; i++) {
132:     if      (i == 0)   flambda[0]   =    h*u[0]   + 2.*d*lambda[0]   - d*lambda[1];
133:     else if (i == 1)   flambda[1]   = 2.*h*u[1]   + 2.*d*lambda[1]   - d*lambda[2];
134:     else if (i == N-1) flambda[N-1] =    h*u[N-1] + 2.*d*lambda[N-1] - d*lambda[N-2];
135:     else if (i == N-2) flambda[N-2] = 2.*h*u[N-2] + 2.*d*lambda[N-2] - d*lambda[N-3];
136:     else               flambda[i]   = 2.*h*u[i]   - d*(lambda[i+1] - 2.0*lambda[i] + lambda[i-1]);
137:   }

139:   /* derivative of L() w.r.t. lambda */
140:   for (i=xs; i<xs+xm; i++) {
141:     if      (i == 0)   fu[0]   = 2.0*d*(u[0] - w[0]);
142:     else if (i == N-1) fu[N-1] = 2.0*d*u[N-1];
143:     else               fu[i]   = -(d*(u[i+1] - 2.0*u[i] + u[i-1]) - 2.0*h);
144:   }

146:   DAVecRestoreArray(user->da1,vu,(void**)&u);
147:   DAVecRestoreArray(user->da1,vfu,(void**)&fu);
148:   DAVecRestoreArray(user->da1,vlambda,(void**)&lambda);
149:   DAVecRestoreArray(user->da1,vflambda,(void**)&flambda);

151:   VecPackGather(user->packer,FU,fw,vfu,vflambda);
152:   VecPackRestoreLocalVectors(user->packer,&w,&vu,&vlambda);
153:   VecPackRestoreLocalVectors(user->packer,&fw,&vfu,&vflambda);
154:   return(0);
155: }

157: int Monitor(SNES snes,int its,PetscReal rnorm,void *dummy)
158: {
159:   UserCtx *user = (UserCtx*)dummy;
160:   int     ierr;
161:   PetscScalar  *w;
162:   Vec     u,lambda,U,F;

165:   SNESGetSolution(snes,&U);
166:   VecPackGetAccess(user->packer,U,&w,&u,&lambda);
167:   VecView(u,user->u_viewer);
168:   VecView(lambda,user->lambda_viewer);
169:   VecPackRestoreAccess(user->packer,U,&w,&u,&lambda);

171:   SNESGetFunction(snes,&F,0,0);
172:   VecPackGetAccess(user->packer,F,&w,&u,&lambda);
173:   VecView(u,user->fu_viewer);
174:   VecView(lambda,user->flambda_viewer);
175:   VecPackRestoreAccess(user->packer,F,&w,&u,&lambda);
176:   return(0);
177: }