xref: /petsc/src/tao/constrained/tutorials/maros.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 /*T
12c4762a1bSJed Brown    Concepts: TAO^Solving an unconstrained minimization problem
13c4762a1bSJed Brown    Routines: TaoCreate(); TaoSetType();
14a82e8c82SStefano Zampini    Routines: TaoSetSolution();
15a82e8c82SStefano Zampini    Routines: TaoSetObjectiveAndGradient();
16c4762a1bSJed Brown    Routines: TaoSetEqualityConstraintsRoutine();
17c4762a1bSJed Brown    Routines: TaoSetInequalityConstraintsRoutine();
18c4762a1bSJed Brown    Routines: TaoSetEqualityJacobianRoutine();
19c4762a1bSJed Brown    Routines: TaoSetInequalityJacobianRoutine();
20a82e8c82SStefano Zampini    Routines: TaoSetHessian(); TaoSetFromOptions();
21c4762a1bSJed Brown    Routines: TaoGetKSP(); TaoSolve();
22c4762a1bSJed Brown    Routines: TaoGetConvergedReason();TaoDestroy();
23c4762a1bSJed Brown    Processors: 1
24c4762a1bSJed Brown T*/
25c4762a1bSJed Brown 
26c4762a1bSJed Brown /*
27c4762a1bSJed Brown    User-defined application context - contains data needed by the
28c4762a1bSJed Brown    application-provided call-back routines, FormFunction(),
29c4762a1bSJed Brown    FormGradient(), and FormHessian().
30c4762a1bSJed Brown */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown /*
33c4762a1bSJed Brown    x,d in R^n
34c4762a1bSJed Brown    f in R
35c4762a1bSJed Brown    bin in R^mi
36c4762a1bSJed Brown    beq in R^me
37c4762a1bSJed Brown    Aeq in R^(me x n)
38c4762a1bSJed Brown    Ain in R^(mi x n)
39c4762a1bSJed Brown    H in R^(n x n)
40c4762a1bSJed Brown    min f=(1/2)*x'*H*x + d'*x
41c4762a1bSJed Brown    s.t.  Aeq*x == beq
42c4762a1bSJed Brown          Ain*x >= bin
43c4762a1bSJed Brown */
44c4762a1bSJed Brown typedef struct {
45c4762a1bSJed Brown   char     name[32];
46c4762a1bSJed Brown   PetscInt n; /* Length x */
47c4762a1bSJed Brown   PetscInt me; /* number of equality constraints */
48c4762a1bSJed Brown   PetscInt mi; /* number of inequality constraints */
49c4762a1bSJed Brown   PetscInt m;  /* me+mi */
50c4762a1bSJed Brown   Mat      Aeq,Ain,H;
51c4762a1bSJed Brown   Vec      beq,bin,d;
52c4762a1bSJed Brown } AppCtx;
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /* -------- User-defined Routines --------- */
55c4762a1bSJed Brown 
56c4762a1bSJed Brown PetscErrorCode InitializeProblem(AppCtx*);
57c4762a1bSJed Brown PetscErrorCode DestroyProblem(AppCtx *);
58c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void *);
59c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat, void*);
60c4762a1bSJed Brown PetscErrorCode FormInequalityConstraints(Tao,Vec,Vec,void*);
61c4762a1bSJed Brown PetscErrorCode FormEqualityConstraints(Tao,Vec,Vec,void*);
62c4762a1bSJed Brown PetscErrorCode FormInequalityJacobian(Tao,Vec,Mat,Mat, void*);
63c4762a1bSJed Brown PetscErrorCode FormEqualityJacobian(Tao,Vec,Mat,Mat, void*);
64c4762a1bSJed Brown 
65c4762a1bSJed Brown PetscErrorCode main(int argc,char **argv)
66c4762a1bSJed Brown {
67c4762a1bSJed Brown   PetscMPIInt        size;
68c4762a1bSJed Brown   Vec                x;                   /* solution */
69c4762a1bSJed Brown   KSP                ksp;
70c4762a1bSJed Brown   PC                 pc;
71c4762a1bSJed Brown   Vec                ceq,cin;
72c4762a1bSJed Brown   PetscBool          flg;                 /* A return value when checking for use options */
73c4762a1bSJed Brown   Tao                tao;                 /* Tao solver context */
74c4762a1bSJed Brown   TaoConvergedReason reason;
75c4762a1bSJed Brown   AppCtx             user;                /* application context */
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Initialize TAO,PETSc */
78*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char *)0,help));
795f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
80c4762a1bSJed Brown   /* Specify default parameters for the problem, check for command-line overrides */
815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(user.name,"HS21",sizeof(user.name)));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-cutername",user.name,sizeof(user.name),&flg));
83c4762a1bSJed Brown 
845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n---- MAROS Problem %s -----\n",user.name));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(InitializeProblem(&user));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user.d,&x));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user.beq,&ceq));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user.bin,&cin));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,1.0));
90c4762a1bSJed Brown 
915f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOIPM));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,x));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*)&user));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetEqualityConstraintsRoutine(tao,ceq,FormEqualityConstraints,(void*)&user));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetInequalityConstraintsRoutine(tao,cin,FormInequalityConstraints,(void*)&user));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetInequalityBounds(tao,user.bin,NULL));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetJacobianEqualityRoutine(tao,user.Aeq,user.Aeq,FormEqualityJacobian,(void*)&user));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetJacobianInequalityRoutine(tao,user.Ain,user.Ain,FormInequalityJacobian,(void*)&user));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetKSP(tao,&ksp));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(ksp,&pc));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(pc,PCLU));
104c4762a1bSJed Brown   /*
105c4762a1bSJed Brown       This algorithm produces matrices with zeros along the diagonal therefore we need to use
106c4762a1bSJed Brown     SuperLU which does partial pivoting
107c4762a1bSJed Brown   */
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetType(ksp,KSPPREONLY));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetTolerances(tao,0,0,0));
111c4762a1bSJed Brown 
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetConvergedReason(tao,&reason));
115c4762a1bSJed Brown   if (reason < 0) {
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge due to %s.\n",TaoConvergedReasons[reason]));
117c4762a1bSJed Brown   } else {
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(MPI_COMM_WORLD, "Optimization completed with status %s.\n",TaoConvergedReasons[reason]));
119c4762a1bSJed Brown   }
120c4762a1bSJed Brown 
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(DestroyProblem(&user));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ceq));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&cin));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
126c4762a1bSJed Brown 
127*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
128*b122ec5aSJacob Faibussowitsch   return 0;
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown PetscErrorCode InitializeProblem(AppCtx *user)
132c4762a1bSJed Brown {
133c4762a1bSJed Brown   PetscViewer loader;
134c4762a1bSJed Brown   MPI_Comm    comm;
135c4762a1bSJed Brown   PetscInt    nrows,ncols,i;
136c4762a1bSJed Brown   PetscScalar one = 1.0;
137c4762a1bSJed Brown   char        filebase[128];
138c4762a1bSJed Brown   char        filename[128];
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   PetscFunctionBegin;
141c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(filebase,user->name,sizeof(filebase)));
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(filebase,"/",sizeof(filebase)));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(filename,filebase,sizeof(filename)));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(filename,"f",sizeof(filename)));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
147c4762a1bSJed Brown 
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(comm,&user->d));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecLoad(user->d,loader));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&loader));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(user->d,&nrows));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user->d));
153c4762a1bSJed Brown   user->n = nrows;
154c4762a1bSJed Brown 
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(filename,filebase,sizeof(filename)));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(filename,"H",sizeof(filename)));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
158c4762a1bSJed Brown 
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(comm,&user->H));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->H,PETSC_DECIDE,PETSC_DECIDE,nrows,nrows));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(user->H,loader));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&loader));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(user->H,&nrows,&ncols));
1643c859ba3SBarry Smith   PetscCheck(nrows == user->n,comm,PETSC_ERR_ARG_SIZ,"H: nrows != n");
1653c859ba3SBarry Smith   PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"H: ncols != n");
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->H));
167c4762a1bSJed Brown 
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(filename,filebase,sizeof(filename)));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(filename,"Aeq",sizeof(filename)));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(comm,&user->Aeq));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(user->Aeq,loader));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&loader));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(user->Aeq,&nrows,&ncols));
1753c859ba3SBarry Smith   PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"Aeq ncols != H nrows");
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->Aeq));
177c4762a1bSJed Brown   user->me = nrows;
178c4762a1bSJed Brown 
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(filename,filebase,sizeof(filename)));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(filename,"Beq",sizeof(filename)));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(comm,&user->beq));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecLoad(user->beq,loader));
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&loader));
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(user->beq,&nrows));
1863c859ba3SBarry Smith   PetscCheck(nrows == user->me,comm,PETSC_ERR_ARG_SIZ,"Aeq nrows != Beq n");
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user->beq));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   user->mi = user->n;
190c4762a1bSJed Brown   /* Ain = eye(n,n) */
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(comm,&user->Ain));
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(user->Ain,MATAIJ));
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->Ain,PETSC_DECIDE,PETSC_DECIDE,user->mi,user->mi));
194c4762a1bSJed Brown 
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(user->Ain,1,NULL,0,NULL));
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(user->Ain,1,NULL));
197c4762a1bSJed Brown 
1985f80ce2aSJacob Faibussowitsch   for (i=0;i<user->mi;i++) CHKERRQ(MatSetValues(user->Ain,1,&i,1,&i,&one,INSERT_VALUES));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(user->Ain,MAT_FINAL_ASSEMBLY));
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(user->Ain,MAT_FINAL_ASSEMBLY));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->Ain));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* bin = [0,0 ... 0]' */
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(comm,&user->bin));
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetType(user->bin,VECMPI));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(user->bin,PETSC_DECIDE,user->mi));
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(user->bin,0.0));
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user->bin));
209c4762a1bSJed Brown   user->m = user->me + user->mi;
210c4762a1bSJed Brown   PetscFunctionReturn(0);
211c4762a1bSJed Brown }
212c4762a1bSJed Brown 
213c4762a1bSJed Brown PetscErrorCode DestroyProblem(AppCtx *user)
214c4762a1bSJed Brown {
215c4762a1bSJed Brown   PetscFunctionBegin;
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->H));
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Aeq));
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Ain));
2195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->beq));
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->bin));
2215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->d));
222c4762a1bSJed Brown   PetscFunctionReturn(0);
223c4762a1bSJed Brown }
2245f80ce2aSJacob Faibussowitsch 
225c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx)
226c4762a1bSJed Brown {
227c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
228c4762a1bSJed Brown   PetscScalar    xtHx;
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   PetscFunctionBegin;
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->H,x,g));
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDot(x,g,&xtHx));
2335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDot(x,user->d,f));
234c4762a1bSJed Brown   *f += 0.5*xtHx;
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(g,1.0,user->d));
236c4762a1bSJed Brown   PetscFunctionReturn(0);
237c4762a1bSJed Brown }
238c4762a1bSJed Brown 
239c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
240c4762a1bSJed Brown {
241c4762a1bSJed Brown   PetscFunctionBegin;
242c4762a1bSJed Brown   PetscFunctionReturn(0);
243c4762a1bSJed Brown }
244c4762a1bSJed Brown 
245c4762a1bSJed Brown PetscErrorCode FormInequalityConstraints(Tao tao, Vec x, Vec ci, void *ctx)
246c4762a1bSJed Brown {
247c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   PetscFunctionBegin;
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->Ain,x,ci));
251c4762a1bSJed Brown   PetscFunctionReturn(0);
252c4762a1bSJed Brown }
253c4762a1bSJed Brown 
254c4762a1bSJed Brown PetscErrorCode FormEqualityConstraints(Tao tao, Vec x, Vec ce,void *ctx)
255c4762a1bSJed Brown {
256c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   PetscFunctionBegin;
2595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->Aeq,x,ce));
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(ce,-1.0,user->beq));
261c4762a1bSJed Brown   PetscFunctionReturn(0);
262c4762a1bSJed Brown }
263c4762a1bSJed Brown 
264c4762a1bSJed Brown PetscErrorCode FormInequalityJacobian(Tao tao, Vec x, Mat JI, Mat JIpre,  void *ctx)
265c4762a1bSJed Brown {
266c4762a1bSJed Brown   PetscFunctionBegin;
267c4762a1bSJed Brown   PetscFunctionReturn(0);
268c4762a1bSJed Brown }
269c4762a1bSJed Brown 
270c4762a1bSJed Brown PetscErrorCode FormEqualityJacobian(Tao tao, Vec x, Mat JE, Mat JEpre, void *ctx)
271c4762a1bSJed Brown {
272c4762a1bSJed Brown   PetscFunctionBegin;
273c4762a1bSJed Brown   PetscFunctionReturn(0);
274c4762a1bSJed Brown }
275c4762a1bSJed Brown 
276c4762a1bSJed Brown /*TEST
277c4762a1bSJed Brown 
278c4762a1bSJed Brown    build:
279c4762a1bSJed Brown       requires: !complex
280c4762a1bSJed Brown 
281c4762a1bSJed Brown    test:
282c4762a1bSJed Brown       requires: superlu
283c4762a1bSJed Brown       localrunfiles: HS21
284c4762a1bSJed Brown 
285c4762a1bSJed Brown TEST*/
286