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