Actual source code: toy.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: /* Program usage: mpiexec -n 1 toy[-help] [all TAO options] */

  3: /* ----------------------------------------------------------------------
  4: min f=(x1-x2)^2 + (x2-2)^2 -2*x1-2*x2
  5: s.t.     x1^2 + x2 = 2
  6:       0 <= x1^2 - x2 <= 1
  7:       -1 <= x1,x2 <= 2
  8: ---------------------------------------------------------------------- */

 10:  #include <petsctao.h>

 12: static  char help[]="";

 14: /*
 15:    User-defined application context - contains data needed by the
 16:    application-provided call-back routines, FormFunction(),
 17:    FormGradient(), and FormHessian().
 18: */

 20: /*
 21:    x,d in R^n
 22:    f in R
 23:    bin in R^mi
 24:    beq in R^me
 25:    Aeq in R^(me x n)
 26:    Ain in R^(mi x n)
 27:    H in R^(n x n)
 28:    min f=(1/2)*x'*H*x + d'*x
 29:    s.t.  Aeq*x == beq
 30:          Ain*x >= bin
 31: */
 32: typedef struct {
 33:   PetscInt n; /* Length x */
 34:   PetscInt ne; /* number of equality constraints */
 35:   PetscInt ni; /* number of inequality constraints */
 36:   Vec      x,xl,xu;
 37:   Vec      ce,ci,bl,bu;
 38:   Mat      Ae,Ai,H;
 39: } AppCtx;

 41: /* -------- User-defined Routines --------- */

 43: PetscErrorCode InitializeProblem(AppCtx *);
 44: PetscErrorCode DestroyProblem(AppCtx *);
 45: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void *);
 46: PetscErrorCode FormHessian(Tao,Vec,Mat,Mat, void*);
 47: PetscErrorCode FormInequalityConstraints(Tao,Vec,Vec,void*);
 48: PetscErrorCode FormEqualityConstraints(Tao,Vec,Vec,void*);
 49: PetscErrorCode FormInequalityJacobian(Tao,Vec,Mat,Mat, void*);
 50: PetscErrorCode FormEqualityJacobian(Tao,Vec,Mat,Mat, void*);



 54: PetscErrorCode main(int argc,char **argv)
 55: {
 56:   PetscErrorCode     ierr;                /* used to check for functions returning nonzeros */
 57:   Tao                tao;
 58:   KSP                ksp;
 59:   PC                 pc;
 60:   AppCtx             user;                /* application context */

 62:   PetscInitialize(&argc,&argv,(char *)0,help);if (ierr) return ierr;
 63:   PetscPrintf(PETSC_COMM_WORLD,"\n---- TOY Problem -----\n");
 64:   PetscPrintf(PETSC_COMM_WORLD,"Solution should be f(1,1)=-2\n");
 65:   InitializeProblem(&user);
 66:   TaoCreate(PETSC_COMM_WORLD,&tao);
 67:   TaoSetType(tao,TAOIPM);
 68:   TaoSetInitialVector(tao,user.x);
 69:   TaoSetVariableBounds(tao,user.xl,user.xu);
 70:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*)&user);

 72:   TaoSetEqualityConstraintsRoutine(tao,user.ce,FormEqualityConstraints,(void*)&user);
 73:   TaoSetInequalityConstraintsRoutine(tao,user.ci,FormInequalityConstraints,(void*)&user);

 75:   TaoSetJacobianEqualityRoutine(tao,user.Ae,user.Ae,FormEqualityJacobian,(void*)&user);
 76:   TaoSetJacobianInequalityRoutine(tao,user.Ai,user.Ai,FormInequalityJacobian,(void*)&user);
 77:   TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void*)&user);
 78:   /* TaoSetTolerances(tao,0,0,0); */

 80:   TaoSetFromOptions(tao);

 82:   TaoGetKSP(tao,&ksp);
 83:   KSPGetPC(ksp,&pc);
 84:   PCSetType(pc,PCLU);
 85:   /*
 86:       This algorithm produces matrices with zeros along the diagonal therefore we need to use
 87:     SuperLU which does partial pivoting
 88:   */
 89:   PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU);
 90:   KSPSetType(ksp,KSPPREONLY);
 91:   KSPSetFromOptions(ksp);

 93:   /* TaoSetTolerances(tao,0,0,0); */
 94:   TaoSolve(tao);

 96:   DestroyProblem(&user);
 97:   TaoDestroy(&tao);
 98:   PetscFinalize();
 99:   return ierr;
100: }

102: PetscErrorCode InitializeProblem(AppCtx *user)
103: {

107:   user->n = 2;
108:   VecCreateSeq(PETSC_COMM_SELF,user->n,&user->x);
109:   VecDuplicate(user->x,&user->xl);
110:   VecDuplicate(user->x,&user->xu);
111:   VecSet(user->x,0.0);
112:   VecSet(user->xl,-1.0);
113:   VecSet(user->xu,2.0);

115:   user->ne = 1;
116:   VecCreateSeq(PETSC_COMM_SELF,user->ne,&user->ce);

118:   user->ni = 2;
119:   VecCreateSeq(PETSC_COMM_SELF,user->ni,&user->ci);

121:   MatCreateSeqAIJ(PETSC_COMM_SELF,user->ne,user->n,user->n,NULL,&user->Ae);
122:   MatCreateSeqAIJ(PETSC_COMM_SELF,user->ni,user->n,user->n,NULL,&user->Ai);
123:   MatSetFromOptions(user->Ae);
124:   MatSetFromOptions(user->Ai);


127:   MatCreateSeqAIJ(PETSC_COMM_SELF,user->n,user->n,1,NULL,&user->H);
128:   MatSetFromOptions(user->H);

130:   return(0);
131: }

133: PetscErrorCode DestroyProblem(AppCtx *user)
134: {

138:   MatDestroy(&user->Ae);
139:   MatDestroy(&user->Ai);
140:   MatDestroy(&user->H);

142:   VecDestroy(&user->x);
143:   VecDestroy(&user->ce);
144:   VecDestroy(&user->ci);
145:   VecDestroy(&user->xl);
146:   VecDestroy(&user->xu);
147:   return(0);
148: }

150: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ctx)
151: {
152:   PetscScalar       *g;
153:   const PetscScalar *x;
154:   PetscErrorCode    ierr;

157:   VecGetArrayRead(X,&x);
158:   VecGetArray(G,&g);
159:   *f = (x[0]-2.0)*(x[0]-2.0) + (x[1]-2.0)*(x[1]-2.0) - 2.0*(x[0]+x[1]);
160:   g[0] = 2.0*(x[0]-2.0) - 2.0;
161:   g[1] = 2.0*(x[1]-2.0) - 2.0;
162:   VecRestoreArrayRead(X,&x);
163:   VecRestoreArray(G,&g);
164:   return(0);
165: }

167: PetscErrorCode FormHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
168: {
169:   Vec               DE,DI;
170:   const PetscScalar *de, *di;
171:   PetscInt          zero=0,one=1;
172:   PetscScalar       two=2.0;
173:   PetscScalar       val;
174:   PetscErrorCode    ierr;

177:   TaoGetDualVariables(tao,&DE,&DI);

179:   VecGetArrayRead(DE,&de);
180:   VecGetArrayRead(DI,&di);
181:   val=2.0 * (1 + de[0] + di[0] - di[1]);
182:   VecRestoreArrayRead(DE,&de);
183:   VecRestoreArrayRead(DI,&di);

185:   MatSetValues(H,1,&zero,1,&zero,&val,INSERT_VALUES);
186:   MatSetValues(H,1,&one,1,&one,&two,INSERT_VALUES);

188:   MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
189:   MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
190:   return(0);
191: }

193: PetscErrorCode FormInequalityConstraints(Tao tao, Vec X, Vec CI, void *ctx)
194: {
195:   const PetscScalar *x;
196:   PetscScalar       *c;
197:   PetscErrorCode    ierr;

200:   VecGetArrayRead(X,&x);
201:   VecGetArray(CI,&c);
202:   c[0] = x[0]*x[0] - x[1];
203:   c[1] = -x[0]*x[0] + x[1] + 1.0;
204:   VecRestoreArrayRead(X,&x);
205:   VecRestoreArray(CI,&c);
206:   return(0);
207: }

209: PetscErrorCode FormEqualityConstraints(Tao tao, Vec X, Vec CE,void *ctx)
210: {
211:   PetscScalar    *x,*c;

215:   VecGetArray(X,&x);
216:   VecGetArray(CE,&c);
217:   c[0] = x[0]*x[0] + x[1] - 2.0;
218:   VecRestoreArray(X,&x);
219:   VecRestoreArray(CE,&c);
220:   return(0);
221: }

223: PetscErrorCode FormInequalityJacobian(Tao tao, Vec X, Mat JI, Mat JIpre,  void *ctx)
224: {
225:   PetscInt          rows[2];
226:   PetscInt          cols[2];
227:   PetscScalar       vals[4];
228:   const PetscScalar *x;
229:   PetscErrorCode    ierr;

232:   VecGetArrayRead(X,&x);
233:   rows[0] = 0;       rows[1] = 1;
234:   cols[0] = 0;       cols[1] = 1;
235:   vals[0] = +2*x[0]; vals[1] = -1.0;
236:   vals[2] = -2*x[0]; vals[3] = +1.0;
237:   VecRestoreArrayRead(X,&x);
238:   MatSetValues(JI,2,rows,2,cols,vals,INSERT_VALUES);
239:   MatAssemblyBegin(JI,MAT_FINAL_ASSEMBLY);
240:   MatAssemblyEnd(JI,MAT_FINAL_ASSEMBLY);

242:   return(0);
243: }

245: PetscErrorCode FormEqualityJacobian(Tao tao, Vec X, Mat JE, Mat JEpre, void *ctx)
246: {
247:   PetscInt          rows[2];
248:   PetscScalar       vals[2];
249:   const PetscScalar *x;
250:   PetscErrorCode    ierr;

253:   VecGetArrayRead(X,&x);
254:   rows[0] = 0;       rows[1] = 1;
255:   vals[0] = 2*x[0];  vals[1] = 1.0;
256:   VecRestoreArrayRead(X,&x);
257:   MatSetValues(JE,1,rows,2,rows,vals,INSERT_VALUES);
258:   MatAssemblyBegin(JE,MAT_FINAL_ASSEMBLY);
259:   MatAssemblyEnd(JE,MAT_FINAL_ASSEMBLY);
260:   return(0);
261: }


264: /*TEST

266:    build:
267:       requires: !complex !define(PETSC_USE_CXX)

269:    test:
270:       requires: superlu
271:       args: -tao_smonitor -tao_view -tao_gatol 1.e-5

273: TEST*/