Actual source code: toy.c
petsc-dev 2014-02-02
1: /* Program usage: mpirun -np 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>
13: static char help[]="";
16: /*
17: User-defined application context - contains data needed by the
18: application-provided call-back routines, FormFunction(),
19: FormGradient(), and FormHessian().
20: */
22: /*
23: x,d in R^n
24: f in R
25: bin in R^mi
26: beq in R^me
27: Aeq in R^(me x n)
28: Ain in R^(mi x n)
29: H in R^(n x n)
30: min f=(1/2)*x'*H*x + d'*x
31: s.t. Aeq*x == beq
32: Ain*x >= bin
33: */
34: typedef struct {
35: PetscInt n; /* Length x */
36: PetscInt ne; /* number of equality constraints */
37: PetscInt ni; /* number of inequality constraints */
38: Vec x,xl,xu;
39: Vec ce,ci,bl,bu;
40: Mat Ae,Ai,H;
42: } AppCtx;
44: /* -------- User-defined Routines --------- */
46: PetscErrorCode InitializeProblem(AppCtx *);
47: PetscErrorCode DestroyProblem(AppCtx *);
48: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void *);
49: PetscErrorCode FormHessian(Tao,Vec,Mat*,Mat*, MatStructure *,void*);
50: PetscErrorCode FormInequalityConstraints(Tao,Vec,Vec,void*);
51: PetscErrorCode FormEqualityConstraints(Tao,Vec,Vec,void*);
52: PetscErrorCode FormInequalityJacobian(Tao,Vec,Mat*,Mat*, MatStructure *,void*);
53: PetscErrorCode FormEqualityJacobian(Tao,Vec,Mat*,Mat*, MatStructure *,void*);
59: PetscErrorCode main(int argc,char **argv)
60: {
62: Tao tao;
63: TaoTerminationReason reason;
64: KSP ksp;
65: PC pc;
66: AppCtx user; /* application context */
68: /* Initialize TAO,PETSc */
69: PetscInitialize(&argc,&argv,(char *)0,help);
71: PetscPrintf(PETSC_COMM_WORLD,"\n---- TOY Problem -----\n");
72: PetscPrintf(PETSC_COMM_WORLD,"Solution should be f(1,1)=-2\n");
73: InitializeProblem(&user);
74: TaoCreate(PETSC_COMM_WORLD,&tao);
75: TaoSetType(tao,"tao_ipm");
76: TaoSetInitialVector(tao,user.x);
77: TaoSetVariableBounds(tao,user.xl,user.xu);
78: TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*)&user);
80: TaoSetEqualityConstraintsRoutine(tao,user.ce,FormEqualityConstraints,(void*)&user);
81: TaoSetInequalityConstraintsRoutine(tao,user.ci,FormInequalityConstraints,(void*)&user);
83: TaoSetJacobianEqualityRoutine(tao,user.Ae,user.Ae,FormEqualityJacobian,(void*)&user);
84: TaoSetJacobianInequalityRoutine(tao,user.Ai,user.Ai,FormInequalityJacobian,(void*)&user);
85: TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void*)&user);
86: TaoSetTolerances(tao,1e-12,0,0,0,0);
88: TaoSetFromOptions(tao);
90: TaoGetKSP(tao,&ksp);
91: KSPGetPC(ksp,&pc);
92: PCFactorSetMatSolverPackage(pc,MATSOLVERSUPERLU);
93: /* TODO -- why didn't that work? */
94: PetscOptionsSetValue("-pc_factor_mat_solver_package","superlu");
95: PCSetType(pc,PCLU);
96: KSPSetType(ksp,KSPPREONLY);
97: KSPSetFromOptions(ksp);
99: TaoSetTolerances(tao,1e-12,0,0,0,0);
101: /* Solve */
102: TaoSolve(tao);
105: /* Analyze solution */
106: TaoGetTerminationReason(tao,&reason);
107: if (reason < 0) {
108: PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge.\n");
109: } else {
110: PetscPrintf(MPI_COMM_WORLD, "Optimization terminated with status %D.\n", reason);
111: }
113: /* Finalize Memory */
114: DestroyProblem(&user);
115: TaoDestroy(&tao);
117: PetscFinalize();
118: return 0;
119: }
123: PetscErrorCode InitializeProblem(AppCtx *user)
124: {
128: user->n = 2;
129: VecCreateSeq(PETSC_COMM_SELF,user->n,&user->x);
130: VecDuplicate(user->x,&user->xl);
131: VecDuplicate(user->x,&user->xu);
132: VecSet(user->x,0.0);
133: VecSet(user->xl,-1.0);
134: VecSet(user->xu,2.0);
136: user->ne = 1;
137: VecCreateSeq(PETSC_COMM_SELF,user->ne,&user->ce);
139: user->ni = 2;
140: VecCreateSeq(PETSC_COMM_SELF,user->ni,&user->ci);
142: MatCreateSeqAIJ(PETSC_COMM_SELF,user->ne,user->n,user->n,NULL,&user->Ae);
143: MatCreateSeqAIJ(PETSC_COMM_SELF,user->ni,user->n,user->n,NULL,&user->Ai);
144: MatSetFromOptions(user->Ae);
145: MatSetFromOptions(user->Ai);
148: MatCreateSeqAIJ(PETSC_COMM_SELF,user->n,user->n,1,NULL,&user->H); MatSetFromOptions(user->H);
151: return(0);
152: }
156: PetscErrorCode DestroyProblem(AppCtx *user)
157: {
161: MatDestroy(&user->Ae);
162: MatDestroy(&user->Ai);
163: MatDestroy(&user->H);
165: VecDestroy(&user->x);
166: VecDestroy(&user->ce);
167: VecDestroy(&user->ci);
168: VecDestroy(&user->xl);
169: VecDestroy(&user->xu);
170: return(0);
171: }
175: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ctx)
176: {
177: PetscScalar *x,*g;
181: VecGetArray(X,&x);
182: VecGetArray(G,&g);
183: *f = (x[0]-2.0)*(x[0]-2.0) + (x[1]-2.0)*(x[1]-2.0) - 2.0*(x[0]+x[1]);
184: g[0] = 2.0*(x[0]-2.0) - 2.0;
185: g[1] = 2.0*(x[1]-2.0) - 2.0;
186: VecRestoreArray(X,&x);
187: VecRestoreArray(G,&g);
188: return(0);
189: }
193: PetscErrorCode FormHessian(Tao tao, Vec x, Mat *H, Mat *Hpre, MatStructure *ms, void *ctx)
194: {
195: AppCtx *user = (AppCtx*)ctx;
196: Vec DE,DI;
197: PetscScalar *de, *di;
198: PetscInt zero=0,one=1;
199: PetscScalar two=2.0;
200: PetscScalar val;
203: val = 0.0;
204: TaoGetDualVariables(tao,&DE,&DI);
206: VecGetArray(DE,&de);
207: VecGetArray(DI,&di);
208: val=2.0 * (1 + de[0] + di[0] - di[1]);
209: VecRestoreArray(DE,&de);
210: VecRestoreArray(DI,&di);
212: MatSetValues(*H,1,&zero,1,&zero,&val,INSERT_VALUES);
213: MatSetValues(*H,1,&one,1,&one,&two,INSERT_VALUES);
215: MatAssemblyBegin(user->H,MAT_FINAL_ASSEMBLY);
216: MatAssemblyEnd(user->H,MAT_FINAL_ASSEMBLY);
218: *ms = SAME_NONZERO_PATTERN;
219: return(0);
220: }
224: PetscErrorCode FormInequalityConstraints(Tao tao, Vec X, Vec CI, void *ctx)
225: {
226: PetscScalar *x,*c;
230: VecGetArray(X,&x);
231: VecGetArray(CI,&c);
232: c[0] = x[0]*x[0] - x[1];
233: c[1] = -x[0]*x[0] + x[1] + 1.0;
234: VecRestoreArray(X,&x);
235: VecRestoreArray(CI,&c);
236: return(0);
238: }
242: PetscErrorCode FormEqualityConstraints(Tao tao, Vec X, Vec CE,void *ctx)
243: {
244: PetscScalar *x,*c;
248: VecGetArray(X,&x);
249: VecGetArray(CE,&c);
250: c[0] = x[0]*x[0] + x[1] - 2.0;
251: VecRestoreArray(X,&x);
252: VecRestoreArray(CE,&c);
253: return(0);
255: }
259: PetscErrorCode FormInequalityJacobian(Tao tao, Vec X, Mat *JI, Mat *JIpre, MatStructure *ms, void *ctx)
260: {
261: PetscInt rows[2];
262: PetscInt cols[2];
263: PetscScalar vals[4],*x;
267: VecGetArray(X,&x);
268: rows[0] = 0; rows[1] = 1;
269: cols[0] = 0; cols[1] = 1;
270: vals[0] = +2*x[0]; vals[1] = -1.0;
271: vals[2] = -2*x[0]; vals[3] = +1.0;
272: VecRestoreArray(X,&x);
273: MatSetValues(*JI,2,rows,2,cols,vals,INSERT_VALUES);
274: MatAssemblyBegin(*JI,MAT_FINAL_ASSEMBLY);
275: MatAssemblyEnd(*JI,MAT_FINAL_ASSEMBLY);
277: return(0);
278: }
283: PetscErrorCode FormEqualityJacobian(Tao tao, Vec X, Mat *JE, Mat *JEpre, MatStructure *ms, void *ctx)
284: {
285: PetscInt rows[2];
286: PetscScalar vals[2],*x;
290: VecGetArray(X,&x);
291: rows[0] = 0; rows[1] = 1;
292: vals[0] = 2*x[0]; vals[1] = 1.0;
293: VecRestoreArray(X,&x);
294: MatSetValues(*JE,1,rows,2,rows,vals,INSERT_VALUES);
295: MatAssemblyBegin(*JE,MAT_FINAL_ASSEMBLY);
296: MatAssemblyEnd(*JE,MAT_FINAL_ASSEMBLY);
298: return(0);
299: }