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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char *)0,help)); 79*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 80c4762a1bSJed Brown /* Specify default parameters for the problem, check for command-line overrides */ 81*9566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(user.name,"HS21",sizeof(user.name))); 82*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-cutername",user.name,sizeof(user.name),&flg)); 83c4762a1bSJed Brown 84*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n---- MAROS Problem %s -----\n",user.name)); 85*9566063dSJacob Faibussowitsch PetscCall(InitializeProblem(&user)); 86*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.d,&x)); 87*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.beq,&ceq)); 88*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.bin,&cin)); 89*9566063dSJacob Faibussowitsch PetscCall(VecSet(x,1.0)); 90c4762a1bSJed Brown 91*9566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao)); 92*9566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOIPM)); 93*9566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 94*9566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*)&user)); 95*9566063dSJacob Faibussowitsch PetscCall(TaoSetEqualityConstraintsRoutine(tao,ceq,FormEqualityConstraints,(void*)&user)); 96*9566063dSJacob Faibussowitsch PetscCall(TaoSetInequalityConstraintsRoutine(tao,cin,FormInequalityConstraints,(void*)&user)); 97*9566063dSJacob Faibussowitsch PetscCall(TaoSetInequalityBounds(tao,user.bin,NULL)); 98*9566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianEqualityRoutine(tao,user.Aeq,user.Aeq,FormEqualityJacobian,(void*)&user)); 99*9566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianInequalityRoutine(tao,user.Ain,user.Ain,FormInequalityJacobian,(void*)&user)); 100*9566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user)); 101*9566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao,&ksp)); 102*9566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 103*9566063dSJacob Faibussowitsch PetscCall(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 */ 108*9566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU)); 109*9566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp,KSPPREONLY)); 110*9566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(tao,0,0,0)); 111c4762a1bSJed Brown 112*9566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 113*9566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 114*9566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(tao,&reason)); 115c4762a1bSJed Brown if (reason < 0) { 116*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge due to %s.\n",TaoConvergedReasons[reason])); 117c4762a1bSJed Brown } else { 118*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "Optimization completed with status %s.\n",TaoConvergedReasons[reason])); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121*9566063dSJacob Faibussowitsch PetscCall(DestroyProblem(&user)); 122*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 123*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ceq)); 124*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cin)); 125*9566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 126c4762a1bSJed Brown 127*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 128b122ec5aSJacob 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; 142*9566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(filebase,user->name,sizeof(filebase))); 143*9566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(filebase,"/",sizeof(filebase))); 144*9566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(filename,filebase,sizeof(filename))); 145*9566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(filename,"f",sizeof(filename))); 146*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader)); 147c4762a1bSJed Brown 148*9566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&user->d)); 149*9566063dSJacob Faibussowitsch PetscCall(VecLoad(user->d,loader)); 150*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&loader)); 151*9566063dSJacob Faibussowitsch PetscCall(VecGetSize(user->d,&nrows)); 152*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->d)); 153c4762a1bSJed Brown user->n = nrows; 154c4762a1bSJed Brown 155*9566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(filename,filebase,sizeof(filename))); 156*9566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(filename,"H",sizeof(filename))); 157*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader)); 158c4762a1bSJed Brown 159*9566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&user->H)); 160*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->H,PETSC_DECIDE,PETSC_DECIDE,nrows,nrows)); 161*9566063dSJacob Faibussowitsch PetscCall(MatLoad(user->H,loader)); 162*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&loader)); 163*9566063dSJacob Faibussowitsch PetscCall(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"); 166*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->H)); 167c4762a1bSJed Brown 168*9566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(filename,filebase,sizeof(filename))); 169*9566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(filename,"Aeq",sizeof(filename))); 170*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader)); 171*9566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&user->Aeq)); 172*9566063dSJacob Faibussowitsch PetscCall(MatLoad(user->Aeq,loader)); 173*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&loader)); 174*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(user->Aeq,&nrows,&ncols)); 1753c859ba3SBarry Smith PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"Aeq ncols != H nrows"); 176*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Aeq)); 177c4762a1bSJed Brown user->me = nrows; 178c4762a1bSJed Brown 179*9566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(filename,filebase,sizeof(filename))); 180*9566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(filename,"Beq",sizeof(filename))); 181*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader)); 182*9566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&user->beq)); 183*9566063dSJacob Faibussowitsch PetscCall(VecLoad(user->beq,loader)); 184*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&loader)); 185*9566063dSJacob Faibussowitsch PetscCall(VecGetSize(user->beq,&nrows)); 1863c859ba3SBarry Smith PetscCheck(nrows == user->me,comm,PETSC_ERR_ARG_SIZ,"Aeq nrows != Beq n"); 187*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->beq)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown user->mi = user->n; 190c4762a1bSJed Brown /* Ain = eye(n,n) */ 191*9566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&user->Ain)); 192*9566063dSJacob Faibussowitsch PetscCall(MatSetType(user->Ain,MATAIJ)); 193*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Ain,PETSC_DECIDE,PETSC_DECIDE,user->mi,user->mi)); 194c4762a1bSJed Brown 195*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Ain,1,NULL,0,NULL)); 196*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Ain,1,NULL)); 197c4762a1bSJed Brown 198*9566063dSJacob Faibussowitsch for (i=0;i<user->mi;i++) PetscCall(MatSetValues(user->Ain,1,&i,1,&i,&one,INSERT_VALUES)); 199*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Ain,MAT_FINAL_ASSEMBLY)); 200*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Ain,MAT_FINAL_ASSEMBLY)); 201*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Ain)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* bin = [0,0 ... 0]' */ 204*9566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&user->bin)); 205*9566063dSJacob Faibussowitsch PetscCall(VecSetType(user->bin,VECMPI)); 206*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->bin,PETSC_DECIDE,user->mi)); 207*9566063dSJacob Faibussowitsch PetscCall(VecSet(user->bin,0.0)); 208*9566063dSJacob Faibussowitsch PetscCall(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; 216*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->H)); 217*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Aeq)); 218*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Ain)); 219*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->beq)); 220*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->bin)); 221*9566063dSJacob Faibussowitsch PetscCall(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; 231*9566063dSJacob Faibussowitsch PetscCall(MatMult(user->H,x,g)); 232*9566063dSJacob Faibussowitsch PetscCall(VecDot(x,g,&xtHx)); 233*9566063dSJacob Faibussowitsch PetscCall(VecDot(x,user->d,f)); 234c4762a1bSJed Brown *f += 0.5*xtHx; 235*9566063dSJacob Faibussowitsch PetscCall(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; 250*9566063dSJacob Faibussowitsch PetscCall(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; 259*9566063dSJacob Faibussowitsch PetscCall(MatMult(user->Aeq,x,ce)); 260*9566063dSJacob Faibussowitsch PetscCall(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