xref: /petsc/src/tao/constrained/tutorials/maros.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 maros1 [-help] [all TAO options] */
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /* ----------------------------------------------------------------------
4c4762a1bSJed Brown TODO Explain maros example
5c4762a1bSJed Brown ---------------------------------------------------------------------- */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #include <petsctao.h>
8c4762a1bSJed Brown 
9c4762a1bSJed Brown static  char help[]="";
10c4762a1bSJed Brown 
11c4762a1bSJed Brown /*
12c4762a1bSJed Brown    User-defined application context - contains data needed by the
13c4762a1bSJed Brown    application-provided call-back routines, FormFunction(),
14c4762a1bSJed Brown    FormGradient(), and FormHessian().
15c4762a1bSJed Brown */
16c4762a1bSJed Brown 
17c4762a1bSJed Brown /*
18c4762a1bSJed Brown    x,d in R^n
19c4762a1bSJed Brown    f in R
20c4762a1bSJed Brown    bin in R^mi
21c4762a1bSJed Brown    beq in R^me
22c4762a1bSJed Brown    Aeq in R^(me x n)
23c4762a1bSJed Brown    Ain in R^(mi x n)
24c4762a1bSJed Brown    H in R^(n x n)
25c4762a1bSJed Brown    min f=(1/2)*x'*H*x + d'*x
26c4762a1bSJed Brown    s.t.  Aeq*x == beq
27c4762a1bSJed Brown          Ain*x >= bin
28c4762a1bSJed Brown */
29c4762a1bSJed Brown typedef struct {
30c4762a1bSJed Brown   char     name[32];
31c4762a1bSJed Brown   PetscInt n; /* Length x */
32c4762a1bSJed Brown   PetscInt me; /* number of equality constraints */
33c4762a1bSJed Brown   PetscInt mi; /* number of inequality constraints */
34c4762a1bSJed Brown   PetscInt m;  /* me+mi */
35c4762a1bSJed Brown   Mat      Aeq,Ain,H;
36c4762a1bSJed Brown   Vec      beq,bin,d;
37c4762a1bSJed Brown } AppCtx;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown /* -------- User-defined Routines --------- */
40c4762a1bSJed Brown 
41c4762a1bSJed Brown PetscErrorCode InitializeProblem(AppCtx*);
42c4762a1bSJed Brown PetscErrorCode DestroyProblem(AppCtx *);
43c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void *);
44c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat, void*);
45c4762a1bSJed Brown PetscErrorCode FormInequalityConstraints(Tao,Vec,Vec,void*);
46c4762a1bSJed Brown PetscErrorCode FormEqualityConstraints(Tao,Vec,Vec,void*);
47c4762a1bSJed Brown PetscErrorCode FormInequalityJacobian(Tao,Vec,Mat,Mat, void*);
48c4762a1bSJed Brown PetscErrorCode FormEqualityJacobian(Tao,Vec,Mat,Mat, void*);
49c4762a1bSJed Brown 
50c4762a1bSJed Brown PetscErrorCode main(int argc,char **argv)
51c4762a1bSJed Brown {
52c4762a1bSJed Brown   PetscMPIInt        size;
53c4762a1bSJed Brown   Vec                x;                   /* solution */
54c4762a1bSJed Brown   KSP                ksp;
55c4762a1bSJed Brown   PC                 pc;
56c4762a1bSJed Brown   Vec                ceq,cin;
57c4762a1bSJed Brown   PetscBool          flg;                 /* A return value when checking for use options */
58c4762a1bSJed Brown   Tao                tao;                 /* Tao solver context */
59c4762a1bSJed Brown   TaoConvergedReason reason;
60c4762a1bSJed Brown   AppCtx             user;                /* application context */
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* Initialize TAO,PETSc */
63*327415f7SBarry Smith   PetscFunctionBeginUser;
649566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char *)0,help));
659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
66c4762a1bSJed Brown   /* Specify default parameters for the problem, check for command-line overrides */
679566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(user.name,"HS21",sizeof(user.name)));
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-cutername",user.name,sizeof(user.name),&flg));
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n---- MAROS Problem %s -----\n",user.name));
719566063dSJacob Faibussowitsch   PetscCall(InitializeProblem(&user));
729566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.d,&x));
739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.beq,&ceq));
749566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.bin,&cin));
759566063dSJacob Faibussowitsch   PetscCall(VecSet(x,1.0));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
789566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOIPM));
799566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,x));
809566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*)&user));
819566063dSJacob Faibussowitsch   PetscCall(TaoSetEqualityConstraintsRoutine(tao,ceq,FormEqualityConstraints,(void*)&user));
829566063dSJacob Faibussowitsch   PetscCall(TaoSetInequalityConstraintsRoutine(tao,cin,FormInequalityConstraints,(void*)&user));
839566063dSJacob Faibussowitsch   PetscCall(TaoSetInequalityBounds(tao,user.bin,NULL));
849566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianEqualityRoutine(tao,user.Aeq,user.Aeq,FormEqualityJacobian,(void*)&user));
859566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianInequalityRoutine(tao,user.Ain,user.Ain,FormInequalityJacobian,(void*)&user));
869566063dSJacob Faibussowitsch   PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user));
879566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao,&ksp));
889566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp,&pc));
899566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCLU));
90c4762a1bSJed Brown   /*
91c4762a1bSJed Brown       This algorithm produces matrices with zeros along the diagonal therefore we need to use
92c4762a1bSJed Brown     SuperLU which does partial pivoting
93c4762a1bSJed Brown   */
949566063dSJacob Faibussowitsch   PetscCall(PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU));
959566063dSJacob Faibussowitsch   PetscCall(KSPSetType(ksp,KSPPREONLY));
969566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(tao,0,0,0));
97c4762a1bSJed Brown 
989566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
999566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
1009566063dSJacob Faibussowitsch   PetscCall(TaoGetConvergedReason(tao,&reason));
101c4762a1bSJed Brown   if (reason < 0) {
1029566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge due to %s.\n",TaoConvergedReasons[reason]));
103c4762a1bSJed Brown   } else {
1049566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(MPI_COMM_WORLD, "Optimization completed with status %s.\n",TaoConvergedReasons[reason]));
105c4762a1bSJed Brown   }
106c4762a1bSJed Brown 
1079566063dSJacob Faibussowitsch   PetscCall(DestroyProblem(&user));
1089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ceq));
1109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cin));
1119566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
112c4762a1bSJed Brown 
1139566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
114b122ec5aSJacob Faibussowitsch   return 0;
115c4762a1bSJed Brown }
116c4762a1bSJed Brown 
117c4762a1bSJed Brown PetscErrorCode InitializeProblem(AppCtx *user)
118c4762a1bSJed Brown {
119c4762a1bSJed Brown   PetscViewer loader;
120c4762a1bSJed Brown   MPI_Comm    comm;
121c4762a1bSJed Brown   PetscInt    nrows,ncols,i;
122c4762a1bSJed Brown   PetscScalar one = 1.0;
123c4762a1bSJed Brown   char        filebase[128];
124c4762a1bSJed Brown   char        filename[128];
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   PetscFunctionBegin;
127c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1289566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(filebase,user->name,sizeof(filebase)));
1299566063dSJacob Faibussowitsch   PetscCall(PetscStrlcat(filebase,"/",sizeof(filebase)));
1309566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(filename,filebase,sizeof(filename)));
1319566063dSJacob Faibussowitsch   PetscCall(PetscStrlcat(filename,"f",sizeof(filename)));
1329566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
133c4762a1bSJed Brown 
1349566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm,&user->d));
1359566063dSJacob Faibussowitsch   PetscCall(VecLoad(user->d,loader));
1369566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&loader));
1379566063dSJacob Faibussowitsch   PetscCall(VecGetSize(user->d,&nrows));
1389566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->d));
139c4762a1bSJed Brown   user->n = nrows;
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(filename,filebase,sizeof(filename)));
1429566063dSJacob Faibussowitsch   PetscCall(PetscStrlcat(filename,"H",sizeof(filename)));
1439566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
144c4762a1bSJed Brown 
1459566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,&user->H));
1469566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->H,PETSC_DECIDE,PETSC_DECIDE,nrows,nrows));
1479566063dSJacob Faibussowitsch   PetscCall(MatLoad(user->H,loader));
1489566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&loader));
1499566063dSJacob Faibussowitsch   PetscCall(MatGetSize(user->H,&nrows,&ncols));
1503c859ba3SBarry Smith   PetscCheck(nrows == user->n,comm,PETSC_ERR_ARG_SIZ,"H: nrows != n");
1513c859ba3SBarry Smith   PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"H: ncols != n");
1529566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->H));
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(filename,filebase,sizeof(filename)));
1559566063dSJacob Faibussowitsch   PetscCall(PetscStrlcat(filename,"Aeq",sizeof(filename)));
1569566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
1579566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,&user->Aeq));
1589566063dSJacob Faibussowitsch   PetscCall(MatLoad(user->Aeq,loader));
1599566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&loader));
1609566063dSJacob Faibussowitsch   PetscCall(MatGetSize(user->Aeq,&nrows,&ncols));
1613c859ba3SBarry Smith   PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"Aeq ncols != H nrows");
1629566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Aeq));
163c4762a1bSJed Brown   user->me = nrows;
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(filename,filebase,sizeof(filename)));
1669566063dSJacob Faibussowitsch   PetscCall(PetscStrlcat(filename,"Beq",sizeof(filename)));
1679566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
1689566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm,&user->beq));
1699566063dSJacob Faibussowitsch   PetscCall(VecLoad(user->beq,loader));
1709566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&loader));
1719566063dSJacob Faibussowitsch   PetscCall(VecGetSize(user->beq,&nrows));
1723c859ba3SBarry Smith   PetscCheck(nrows == user->me,comm,PETSC_ERR_ARG_SIZ,"Aeq nrows != Beq n");
1739566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->beq));
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   user->mi = user->n;
176c4762a1bSJed Brown   /* Ain = eye(n,n) */
1779566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,&user->Ain));
1789566063dSJacob Faibussowitsch   PetscCall(MatSetType(user->Ain,MATAIJ));
1799566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Ain,PETSC_DECIDE,PETSC_DECIDE,user->mi,user->mi));
180c4762a1bSJed Brown 
1819566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Ain,1,NULL,0,NULL));
1829566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Ain,1,NULL));
183c4762a1bSJed Brown 
1849566063dSJacob Faibussowitsch   for (i=0;i<user->mi;i++) PetscCall(MatSetValues(user->Ain,1,&i,1,&i,&one,INSERT_VALUES));
1859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Ain,MAT_FINAL_ASSEMBLY));
1869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Ain,MAT_FINAL_ASSEMBLY));
1879566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Ain));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* bin = [0,0 ... 0]' */
1909566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm,&user->bin));
1919566063dSJacob Faibussowitsch   PetscCall(VecSetType(user->bin,VECMPI));
1929566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->bin,PETSC_DECIDE,user->mi));
1939566063dSJacob Faibussowitsch   PetscCall(VecSet(user->bin,0.0));
1949566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->bin));
195c4762a1bSJed Brown   user->m = user->me + user->mi;
196c4762a1bSJed Brown   PetscFunctionReturn(0);
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown PetscErrorCode DestroyProblem(AppCtx *user)
200c4762a1bSJed Brown {
201c4762a1bSJed Brown   PetscFunctionBegin;
2029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->H));
2039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Aeq));
2049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Ain));
2059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->beq));
2069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->bin));
2079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->d));
208c4762a1bSJed Brown   PetscFunctionReturn(0);
209c4762a1bSJed Brown }
2105f80ce2aSJacob Faibussowitsch 
211c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx)
212c4762a1bSJed Brown {
213c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
214c4762a1bSJed Brown   PetscScalar    xtHx;
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   PetscFunctionBegin;
2179566063dSJacob Faibussowitsch   PetscCall(MatMult(user->H,x,g));
2189566063dSJacob Faibussowitsch   PetscCall(VecDot(x,g,&xtHx));
2199566063dSJacob Faibussowitsch   PetscCall(VecDot(x,user->d,f));
220c4762a1bSJed Brown   *f += 0.5*xtHx;
2219566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g,1.0,user->d));
222c4762a1bSJed Brown   PetscFunctionReturn(0);
223c4762a1bSJed Brown }
224c4762a1bSJed Brown 
225c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
226c4762a1bSJed Brown {
227c4762a1bSJed Brown   PetscFunctionBegin;
228c4762a1bSJed Brown   PetscFunctionReturn(0);
229c4762a1bSJed Brown }
230c4762a1bSJed Brown 
231c4762a1bSJed Brown PetscErrorCode FormInequalityConstraints(Tao tao, Vec x, Vec ci, void *ctx)
232c4762a1bSJed Brown {
233c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   PetscFunctionBegin;
2369566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Ain,x,ci));
237c4762a1bSJed Brown   PetscFunctionReturn(0);
238c4762a1bSJed Brown }
239c4762a1bSJed Brown 
240c4762a1bSJed Brown PetscErrorCode FormEqualityConstraints(Tao tao, Vec x, Vec ce,void *ctx)
241c4762a1bSJed Brown {
242c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   PetscFunctionBegin;
2459566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Aeq,x,ce));
2469566063dSJacob Faibussowitsch   PetscCall(VecAXPY(ce,-1.0,user->beq));
247c4762a1bSJed Brown   PetscFunctionReturn(0);
248c4762a1bSJed Brown }
249c4762a1bSJed Brown 
250c4762a1bSJed Brown PetscErrorCode FormInequalityJacobian(Tao tao, Vec x, Mat JI, Mat JIpre,  void *ctx)
251c4762a1bSJed Brown {
252c4762a1bSJed Brown   PetscFunctionBegin;
253c4762a1bSJed Brown   PetscFunctionReturn(0);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown PetscErrorCode FormEqualityJacobian(Tao tao, Vec x, Mat JE, Mat JEpre, void *ctx)
257c4762a1bSJed Brown {
258c4762a1bSJed Brown   PetscFunctionBegin;
259c4762a1bSJed Brown   PetscFunctionReturn(0);
260c4762a1bSJed Brown }
261c4762a1bSJed Brown 
262c4762a1bSJed Brown /*TEST
263c4762a1bSJed Brown 
264c4762a1bSJed Brown    build:
265c4762a1bSJed Brown       requires: !complex
266c4762a1bSJed Brown 
267c4762a1bSJed Brown    test:
268c4762a1bSJed Brown       requires: superlu
269c4762a1bSJed Brown       localrunfiles: HS21
270c4762a1bSJed Brown 
271c4762a1bSJed Brown TEST*/
272