xref: /petsc/src/tao/constrained/tutorials/maros.c (revision 2e3d3ef9bcea2ec3187269b59709989dccbaffee)
1 /* Program usage: mpiexec -n 1 maros1 [-help] [all TAO options] */
2 
3 /* ----------------------------------------------------------------------
4 TODO Explain maros example
5 ---------------------------------------------------------------------- */
6 
7 #include <petsctao.h>
8 
9 static  char help[]="";
10 
11 /*
12    User-defined application context - contains data needed by the
13    application-provided call-back routines, FormFunction(),
14    FormGradient(), and FormHessian().
15 */
16 
17 /*
18    x,d in R^n
19    f in R
20    bin in R^mi
21    beq in R^me
22    Aeq in R^(me x n)
23    Ain in R^(mi x n)
24    H in R^(n x n)
25    min f=(1/2)*x'*H*x + d'*x
26    s.t.  Aeq*x == beq
27          Ain*x >= bin
28 */
29 typedef struct {
30   char     name[32];
31   PetscInt n; /* Length x */
32   PetscInt me; /* number of equality constraints */
33   PetscInt mi; /* number of inequality constraints */
34   PetscInt m;  /* me+mi */
35   Mat      Aeq,Ain,H;
36   Vec      beq,bin,d;
37 } AppCtx;
38 
39 /* -------- User-defined Routines --------- */
40 
41 PetscErrorCode InitializeProblem(AppCtx*);
42 PetscErrorCode DestroyProblem(AppCtx *);
43 PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void *);
44 PetscErrorCode FormHessian(Tao,Vec,Mat,Mat, void*);
45 PetscErrorCode FormInequalityConstraints(Tao,Vec,Vec,void*);
46 PetscErrorCode FormEqualityConstraints(Tao,Vec,Vec,void*);
47 PetscErrorCode FormInequalityJacobian(Tao,Vec,Mat,Mat, void*);
48 PetscErrorCode FormEqualityJacobian(Tao,Vec,Mat,Mat, void*);
49 
50 PetscErrorCode main(int argc,char **argv)
51 {
52   PetscMPIInt        size;
53   Vec                x;                   /* solution */
54   KSP                ksp;
55   PC                 pc;
56   Vec                ceq,cin;
57   PetscBool          flg;                 /* A return value when checking for use options */
58   Tao                tao;                 /* Tao solver context */
59   TaoConvergedReason reason;
60   AppCtx             user;                /* application context */
61 
62   /* Initialize TAO,PETSc */
63   PetscCall(PetscInitialize(&argc,&argv,(char *)0,help));
64   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
65   /* Specify default parameters for the problem, check for command-line overrides */
66   PetscCall(PetscStrncpy(user.name,"HS21",sizeof(user.name)));
67   PetscCall(PetscOptionsGetString(NULL,NULL,"-cutername",user.name,sizeof(user.name),&flg));
68 
69   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n---- MAROS Problem %s -----\n",user.name));
70   PetscCall(InitializeProblem(&user));
71   PetscCall(VecDuplicate(user.d,&x));
72   PetscCall(VecDuplicate(user.beq,&ceq));
73   PetscCall(VecDuplicate(user.bin,&cin));
74   PetscCall(VecSet(x,1.0));
75 
76   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
77   PetscCall(TaoSetType(tao,TAOIPM));
78   PetscCall(TaoSetSolution(tao,x));
79   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*)&user));
80   PetscCall(TaoSetEqualityConstraintsRoutine(tao,ceq,FormEqualityConstraints,(void*)&user));
81   PetscCall(TaoSetInequalityConstraintsRoutine(tao,cin,FormInequalityConstraints,(void*)&user));
82   PetscCall(TaoSetInequalityBounds(tao,user.bin,NULL));
83   PetscCall(TaoSetJacobianEqualityRoutine(tao,user.Aeq,user.Aeq,FormEqualityJacobian,(void*)&user));
84   PetscCall(TaoSetJacobianInequalityRoutine(tao,user.Ain,user.Ain,FormInequalityJacobian,(void*)&user));
85   PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user));
86   PetscCall(TaoGetKSP(tao,&ksp));
87   PetscCall(KSPGetPC(ksp,&pc));
88   PetscCall(PCSetType(pc,PCLU));
89   /*
90       This algorithm produces matrices with zeros along the diagonal therefore we need to use
91     SuperLU which does partial pivoting
92   */
93   PetscCall(PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU));
94   PetscCall(KSPSetType(ksp,KSPPREONLY));
95   PetscCall(TaoSetTolerances(tao,0,0,0));
96 
97   PetscCall(TaoSetFromOptions(tao));
98   PetscCall(TaoSolve(tao));
99   PetscCall(TaoGetConvergedReason(tao,&reason));
100   if (reason < 0) {
101     PetscCall(PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge due to %s.\n",TaoConvergedReasons[reason]));
102   } else {
103     PetscCall(PetscPrintf(MPI_COMM_WORLD, "Optimization completed with status %s.\n",TaoConvergedReasons[reason]));
104   }
105 
106   PetscCall(DestroyProblem(&user));
107   PetscCall(VecDestroy(&x));
108   PetscCall(VecDestroy(&ceq));
109   PetscCall(VecDestroy(&cin));
110   PetscCall(TaoDestroy(&tao));
111 
112   PetscCall(PetscFinalize());
113   return 0;
114 }
115 
116 PetscErrorCode InitializeProblem(AppCtx *user)
117 {
118   PetscViewer loader;
119   MPI_Comm    comm;
120   PetscInt    nrows,ncols,i;
121   PetscScalar one = 1.0;
122   char        filebase[128];
123   char        filename[128];
124 
125   PetscFunctionBegin;
126   comm = PETSC_COMM_WORLD;
127   PetscCall(PetscStrncpy(filebase,user->name,sizeof(filebase)));
128   PetscCall(PetscStrlcat(filebase,"/",sizeof(filebase)));
129   PetscCall(PetscStrncpy(filename,filebase,sizeof(filename)));
130   PetscCall(PetscStrlcat(filename,"f",sizeof(filename)));
131   PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
132 
133   PetscCall(VecCreate(comm,&user->d));
134   PetscCall(VecLoad(user->d,loader));
135   PetscCall(PetscViewerDestroy(&loader));
136   PetscCall(VecGetSize(user->d,&nrows));
137   PetscCall(VecSetFromOptions(user->d));
138   user->n = nrows;
139 
140   PetscCall(PetscStrncpy(filename,filebase,sizeof(filename)));
141   PetscCall(PetscStrlcat(filename,"H",sizeof(filename)));
142   PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
143 
144   PetscCall(MatCreate(comm,&user->H));
145   PetscCall(MatSetSizes(user->H,PETSC_DECIDE,PETSC_DECIDE,nrows,nrows));
146   PetscCall(MatLoad(user->H,loader));
147   PetscCall(PetscViewerDestroy(&loader));
148   PetscCall(MatGetSize(user->H,&nrows,&ncols));
149   PetscCheck(nrows == user->n,comm,PETSC_ERR_ARG_SIZ,"H: nrows != n");
150   PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"H: ncols != n");
151   PetscCall(MatSetFromOptions(user->H));
152 
153   PetscCall(PetscStrncpy(filename,filebase,sizeof(filename)));
154   PetscCall(PetscStrlcat(filename,"Aeq",sizeof(filename)));
155   PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
156   PetscCall(MatCreate(comm,&user->Aeq));
157   PetscCall(MatLoad(user->Aeq,loader));
158   PetscCall(PetscViewerDestroy(&loader));
159   PetscCall(MatGetSize(user->Aeq,&nrows,&ncols));
160   PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"Aeq ncols != H nrows");
161   PetscCall(MatSetFromOptions(user->Aeq));
162   user->me = nrows;
163 
164   PetscCall(PetscStrncpy(filename,filebase,sizeof(filename)));
165   PetscCall(PetscStrlcat(filename,"Beq",sizeof(filename)));
166   PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
167   PetscCall(VecCreate(comm,&user->beq));
168   PetscCall(VecLoad(user->beq,loader));
169   PetscCall(PetscViewerDestroy(&loader));
170   PetscCall(VecGetSize(user->beq,&nrows));
171   PetscCheck(nrows == user->me,comm,PETSC_ERR_ARG_SIZ,"Aeq nrows != Beq n");
172   PetscCall(VecSetFromOptions(user->beq));
173 
174   user->mi = user->n;
175   /* Ain = eye(n,n) */
176   PetscCall(MatCreate(comm,&user->Ain));
177   PetscCall(MatSetType(user->Ain,MATAIJ));
178   PetscCall(MatSetSizes(user->Ain,PETSC_DECIDE,PETSC_DECIDE,user->mi,user->mi));
179 
180   PetscCall(MatMPIAIJSetPreallocation(user->Ain,1,NULL,0,NULL));
181   PetscCall(MatSeqAIJSetPreallocation(user->Ain,1,NULL));
182 
183   for (i=0;i<user->mi;i++) PetscCall(MatSetValues(user->Ain,1,&i,1,&i,&one,INSERT_VALUES));
184   PetscCall(MatAssemblyBegin(user->Ain,MAT_FINAL_ASSEMBLY));
185   PetscCall(MatAssemblyEnd(user->Ain,MAT_FINAL_ASSEMBLY));
186   PetscCall(MatSetFromOptions(user->Ain));
187 
188   /* bin = [0,0 ... 0]' */
189   PetscCall(VecCreate(comm,&user->bin));
190   PetscCall(VecSetType(user->bin,VECMPI));
191   PetscCall(VecSetSizes(user->bin,PETSC_DECIDE,user->mi));
192   PetscCall(VecSet(user->bin,0.0));
193   PetscCall(VecSetFromOptions(user->bin));
194   user->m = user->me + user->mi;
195   PetscFunctionReturn(0);
196 }
197 
198 PetscErrorCode DestroyProblem(AppCtx *user)
199 {
200   PetscFunctionBegin;
201   PetscCall(MatDestroy(&user->H));
202   PetscCall(MatDestroy(&user->Aeq));
203   PetscCall(MatDestroy(&user->Ain));
204   PetscCall(VecDestroy(&user->beq));
205   PetscCall(VecDestroy(&user->bin));
206   PetscCall(VecDestroy(&user->d));
207   PetscFunctionReturn(0);
208 }
209 
210 PetscErrorCode FormFunctionGradient(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx)
211 {
212   AppCtx         *user = (AppCtx*)ctx;
213   PetscScalar    xtHx;
214 
215   PetscFunctionBegin;
216   PetscCall(MatMult(user->H,x,g));
217   PetscCall(VecDot(x,g,&xtHx));
218   PetscCall(VecDot(x,user->d,f));
219   *f += 0.5*xtHx;
220   PetscCall(VecAXPY(g,1.0,user->d));
221   PetscFunctionReturn(0);
222 }
223 
224 PetscErrorCode FormHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
225 {
226   PetscFunctionBegin;
227   PetscFunctionReturn(0);
228 }
229 
230 PetscErrorCode FormInequalityConstraints(Tao tao, Vec x, Vec ci, void *ctx)
231 {
232   AppCtx         *user = (AppCtx*)ctx;
233 
234   PetscFunctionBegin;
235   PetscCall(MatMult(user->Ain,x,ci));
236   PetscFunctionReturn(0);
237 }
238 
239 PetscErrorCode FormEqualityConstraints(Tao tao, Vec x, Vec ce,void *ctx)
240 {
241   AppCtx         *user = (AppCtx*)ctx;
242 
243   PetscFunctionBegin;
244   PetscCall(MatMult(user->Aeq,x,ce));
245   PetscCall(VecAXPY(ce,-1.0,user->beq));
246   PetscFunctionReturn(0);
247 }
248 
249 PetscErrorCode FormInequalityJacobian(Tao tao, Vec x, Mat JI, Mat JIpre,  void *ctx)
250 {
251   PetscFunctionBegin;
252   PetscFunctionReturn(0);
253 }
254 
255 PetscErrorCode FormEqualityJacobian(Tao tao, Vec x, Mat JE, Mat JEpre, void *ctx)
256 {
257   PetscFunctionBegin;
258   PetscFunctionReturn(0);
259 }
260 
261 /*TEST
262 
263    build:
264       requires: !complex
265 
266    test:
267       requires: superlu
268       localrunfiles: HS21
269 
270 TEST*/
271