Actual source code: toy.c
petsc-3.7.7 2017-09-25
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*);
56: PetscErrorCode main(int argc,char **argv)
57: {
58: PetscErrorCode ierr; /* used to check for functions returning nonzeros */
59: Tao tao;
60: KSP ksp;
61: PC pc;
62: AppCtx user; /* application context */
64: PetscInitialize(&argc,&argv,(char *)0,help);
65: PetscPrintf(PETSC_COMM_WORLD,"\n---- TOY Problem -----\n");
66: PetscPrintf(PETSC_COMM_WORLD,"Solution should be f(1,1)=-2\n");
67: InitializeProblem(&user);
68: TaoCreate(PETSC_COMM_WORLD,&tao);
69: TaoSetType(tao,TAOIPM);
70: TaoSetInitialVector(tao,user.x);
71: TaoSetVariableBounds(tao,user.xl,user.xu);
72: TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*)&user);
74: TaoSetEqualityConstraintsRoutine(tao,user.ce,FormEqualityConstraints,(void*)&user);
75: TaoSetInequalityConstraintsRoutine(tao,user.ci,FormInequalityConstraints,(void*)&user);
77: TaoSetJacobianEqualityRoutine(tao,user.Ae,user.Ae,FormEqualityJacobian,(void*)&user);
78: TaoSetJacobianInequalityRoutine(tao,user.Ai,user.Ai,FormInequalityJacobian,(void*)&user);
79: TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void*)&user);
80: TaoSetTolerances(tao,0,0,0);
82: TaoSetFromOptions(tao);
84: TaoGetKSP(tao,&ksp);
85: KSPGetPC(ksp,&pc);
86: PCSetType(pc,PCLU);
87: /*
88: This algorithm produces matrices with zeros along the diagonal therefore we need to use
89: SuperLU which does partial pivoting
90: */
91: PCFactorSetMatSolverPackage(pc,MATSOLVERSUPERLU);
92: KSPSetType(ksp,KSPPREONLY);
93: KSPSetFromOptions(ksp);
95: TaoSetTolerances(tao,0,0,0);
96: TaoSolve(tao);
98: DestroyProblem(&user);
99: TaoDestroy(&tao);
100: PetscFinalize();
101: return 0;
102: }
106: PetscErrorCode InitializeProblem(AppCtx *user)
107: {
111: user->n = 2;
112: VecCreateSeq(PETSC_COMM_SELF,user->n,&user->x);
113: VecDuplicate(user->x,&user->xl);
114: VecDuplicate(user->x,&user->xu);
115: VecSet(user->x,0.0);
116: VecSet(user->xl,-1.0);
117: VecSet(user->xu,2.0);
119: user->ne = 1;
120: VecCreateSeq(PETSC_COMM_SELF,user->ne,&user->ce);
122: user->ni = 2;
123: VecCreateSeq(PETSC_COMM_SELF,user->ni,&user->ci);
125: MatCreateSeqAIJ(PETSC_COMM_SELF,user->ne,user->n,user->n,NULL,&user->Ae);
126: MatCreateSeqAIJ(PETSC_COMM_SELF,user->ni,user->n,user->n,NULL,&user->Ai);
127: MatSetFromOptions(user->Ae);
128: MatSetFromOptions(user->Ai);
131: MatCreateSeqAIJ(PETSC_COMM_SELF,user->n,user->n,1,NULL,&user->H);
132: MatSetFromOptions(user->H);
134: return(0);
135: }
139: PetscErrorCode DestroyProblem(AppCtx *user)
140: {
144: MatDestroy(&user->Ae);
145: MatDestroy(&user->Ai);
146: MatDestroy(&user->H);
148: VecDestroy(&user->x);
149: VecDestroy(&user->ce);
150: VecDestroy(&user->ci);
151: VecDestroy(&user->xl);
152: VecDestroy(&user->xu);
153: return(0);
154: }
158: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ctx)
159: {
160: PetscScalar *g;
161: const PetscScalar *x;
162: PetscErrorCode ierr;
165: VecGetArrayRead(X,&x);
166: VecGetArray(G,&g);
167: *f = (x[0]-2.0)*(x[0]-2.0) + (x[1]-2.0)*(x[1]-2.0) - 2.0*(x[0]+x[1]);
168: g[0] = 2.0*(x[0]-2.0) - 2.0;
169: g[1] = 2.0*(x[1]-2.0) - 2.0;
170: VecRestoreArrayRead(X,&x);
171: VecRestoreArray(G,&g);
172: return(0);
173: }
177: PetscErrorCode FormHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
178: {
179: Vec DE,DI;
180: const PetscScalar *de, *di;
181: PetscInt zero=0,one=1;
182: PetscScalar two=2.0;
183: PetscScalar val;
184: PetscErrorCode ierr;
187: TaoGetDualVariables(tao,&DE,&DI);
189: VecGetArrayRead(DE,&de);
190: VecGetArrayRead(DI,&di);
191: val=2.0 * (1 + de[0] + di[0] - di[1]);
192: VecRestoreArrayRead(DE,&de);
193: VecRestoreArrayRead(DI,&di);
195: MatSetValues(H,1,&zero,1,&zero,&val,INSERT_VALUES);
196: MatSetValues(H,1,&one,1,&one,&two,INSERT_VALUES);
198: MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);
199: MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);
200: return(0);
201: }
205: PetscErrorCode FormInequalityConstraints(Tao tao, Vec X, Vec CI, void *ctx)
206: {
207: const PetscScalar *x;
208: PetscScalar *c;
209: PetscErrorCode ierr;
212: VecGetArrayRead(X,&x);
213: VecGetArray(CI,&c);
214: c[0] = x[0]*x[0] - x[1];
215: c[1] = -x[0]*x[0] + x[1] + 1.0;
216: VecRestoreArrayRead(X,&x);
217: VecRestoreArray(CI,&c);
218: return(0);
219: }
223: PetscErrorCode FormEqualityConstraints(Tao tao, Vec X, Vec CE,void *ctx)
224: {
225: PetscScalar *x,*c;
229: VecGetArray(X,&x);
230: VecGetArray(CE,&c);
231: c[0] = x[0]*x[0] + x[1] - 2.0;
232: VecRestoreArray(X,&x);
233: VecRestoreArray(CE,&c);
234: return(0);
235: }
239: PetscErrorCode FormInequalityJacobian(Tao tao, Vec X, Mat JI, Mat JIpre, void *ctx)
240: {
241: PetscInt rows[2];
242: PetscInt cols[2];
243: PetscScalar vals[4];
244: const PetscScalar *x;
245: PetscErrorCode ierr;
248: VecGetArrayRead(X,&x);
249: rows[0] = 0; rows[1] = 1;
250: cols[0] = 0; cols[1] = 1;
251: vals[0] = +2*x[0]; vals[1] = -1.0;
252: vals[2] = -2*x[0]; vals[3] = +1.0;
253: VecRestoreArrayRead(X,&x);
254: MatSetValues(JI,2,rows,2,cols,vals,INSERT_VALUES);
255: MatAssemblyBegin(JI,MAT_FINAL_ASSEMBLY);
256: MatAssemblyEnd(JI,MAT_FINAL_ASSEMBLY);
258: return(0);
259: }
263: PetscErrorCode FormEqualityJacobian(Tao tao, Vec X, Mat JE, Mat JEpre, void *ctx)
264: {
265: PetscInt rows[2];
266: PetscScalar vals[2];
267: const PetscScalar *x;
268: PetscErrorCode ierr;
271: VecGetArrayRead(X,&x);
272: rows[0] = 0; rows[1] = 1;
273: vals[0] = 2*x[0]; vals[1] = 1.0;
274: VecRestoreArrayRead(X,&x);
275: MatSetValues(JE,1,rows,2,rows,vals,INSERT_VALUES);
276: MatAssemblyBegin(JE,MAT_FINAL_ASSEMBLY);
277: MatAssemblyEnd(JE,MAT_FINAL_ASSEMBLY);
278: return(0);
279: }