1 #include "taolinesearch.h" 2 #include "ipm.h" /*I "ipm.h" I*/ 3 /* 4 #define DEBUG_IPM 5 #define DEBUG_K 6 #define DEBUG_SCATTER 7 #define DEBUG_KKT 8 */ 9 /* 10 x,d in R^n 11 f in R 12 nb = mi + nlb+nub 13 s in R^nb is slack vector CI(x) / x-XL / -x+XU 14 bin in R^mi (tao->constraints_inequality) 15 beq in R^me (tao->constraints_equality) 16 lamdai in R^nb (ipmP->lamdai) 17 lamdae in R^me (ipmP->lamdae) 18 Jeq in R^(me x n) (tao->jacobian_equality) 19 Jin in R^(mi x n) (tao->jacobian_inequality) 20 Ai in R^(nb x n) (ipmP->Ai) 21 H in R^(n x n) (tao->hessian) 22 min f=(1/2)*x'*H*x + d'*x 23 s.t. CE(x) == 0 24 CI(x) >= 0 25 x >= tao->XL 26 -x >= -tao->XU 27 */ 28 29 static PetscErrorCode IPMComputeKKT(TaoSolver tao); 30 static PetscErrorCode IPMPushInitialPoint(TaoSolver tao); 31 static PetscErrorCode IPMEvaluate(TaoSolver tao); 32 static PetscErrorCode IPMUpdateK(TaoSolver tao); 33 static PetscErrorCode IPMUpdateAi(TaoSolver tao); 34 static PetscErrorCode IPMGatherRHS(TaoSolver tao,Vec,Vec,Vec,Vec,Vec); 35 static PetscErrorCode IPMScatterStep(TaoSolver tao,Vec,Vec,Vec,Vec,Vec); 36 static PetscErrorCode IPMInitializeBounds(TaoSolver tao); 37 38 #undef __FUNCT__ 39 #define __FUNCT__ "TaoSolve_IPM" 40 static PetscErrorCode TaoSolve_IPM(TaoSolver tao) 41 { 42 PetscErrorCode ierr; 43 TAO_IPM* ipmP = (TAO_IPM*)tao->data; 44 45 46 TaoSolverTerminationReason reason = TAO_CONTINUE_ITERATING; 47 PetscInt iter = 0,its,i; 48 PetscScalar stepsize=1.0; 49 PetscScalar step_s,step_l,alpha,tau,sigma,phi_target; 50 PetscFunctionBegin; 51 52 /* Push initial point away from bounds */ 53 #if defined(DEBUG_IPM) 54 ierr = VecNorm(tao->solution,NORM_2,&tau); CHKERRQ(ierr); 55 VecView(tao->solution,0); 56 #endif 57 ierr = IPMInitializeBounds(tao); CHKERRQ(ierr); 58 ierr = IPMPushInitialPoint(tao); CHKERRQ(ierr); 59 #if defined(DEBUG_IPM) 60 ierr = VecNorm(tao->solution,NORM_2,&tau); CHKERRQ(ierr); 61 VecView(tao->solution,0); 62 /* PetscPrintf(PETSC_COMM_WORLD,"||x0|| = %g\n",tau); */ 63 #endif 64 ierr = VecCopy(tao->solution,ipmP->rhs_x); CHKERRQ(ierr); 65 ierr = IPMEvaluate(tao); CHKERRQ(ierr); 66 ierr = IPMComputeKKT(tao); CHKERRQ(ierr); 67 ierr = TaoMonitor(tao,iter++,ipmP->kkt_f,ipmP->phi,0.0,1.0,&reason); 68 69 70 while (reason == TAO_CONTINUE_ITERATING) { 71 ierr = IPMUpdateK(tao); CHKERRQ(ierr); 72 /* 73 rhs.x = -rd 74 rhs.lame = -rpe 75 rhs.lami = -rpi 76 rhs.com = -com 77 */ 78 79 ierr = VecCopy(ipmP->rd,ipmP->rhs_x); CHKERRQ(ierr); 80 if (ipmP->me > 0) { 81 ierr = VecCopy(ipmP->rpe,ipmP->rhs_lamdae); CHKERRQ(ierr); 82 } 83 if (ipmP->nb > 0) { 84 ierr = VecCopy(ipmP->rpi,ipmP->rhs_lamdai); CHKERRQ(ierr); 85 ierr = VecCopy(ipmP->complementarity,ipmP->rhs_s); CHKERRQ(ierr); 86 } 87 ierr = IPMGatherRHS(tao,ipmP->bigrhs,ipmP->rhs_x,ipmP->rhs_lamdae, 88 ipmP->rhs_lamdai,ipmP->rhs_s); CHKERRQ(ierr); 89 ierr = VecScale(ipmP->bigrhs,-1.0); CHKERRQ(ierr); 90 91 /* solve K * step = rhs */ 92 ierr = KSPSetOperators(tao->ksp,ipmP->K,ipmP->K,DIFFERENT_NONZERO_PATTERN); CHKERRQ(ierr); 93 ierr = KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep);CHKERRQ(ierr); 94 95 ierr = IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds, 96 ipmP->dlamdae,ipmP->dlamdai); CHKERRQ(ierr); 97 ierr = KSPGetIterationNumber(tao->ksp,&its); CHKERRQ(ierr); 98 tao->ksp_its += its; 99 #if defined DEBUG_KKT 100 PetscPrintf(PETSC_COMM_WORLD,"first solve.\n"); 101 PetscPrintf(PETSC_COMM_WORLD,"rhs_lamdai\n"); 102 /* VecView(ipmP->rhs_lamdai,0); 103 ierr = VecView(ipmP->bigrhs,0); 104 ierr = VecView(ipmP->bigstep,0); */ 105 PetscScalar norm1,norm2; 106 ierr = VecNorm(ipmP->bigrhs,NORM_2,&norm1); 107 ierr = VecNorm(ipmP->bigstep,NORM_2,&norm2); 108 PetscPrintf(PETSC_COMM_WORLD,"||rhs|| = %g\t ||step|| = %g\n",norm1,norm2); 109 #endif 110 /* Find distance along step direction to closest bound */ 111 if (ipmP->nb > 0) { 112 ierr = VecStepBoundInfo(ipmP->s,ipmP->Zero_nb,ipmP->Inf_nb,ipmP->ds,&step_s,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 113 ierr = VecStepBoundInfo(ipmP->lamdai,ipmP->Zero_nb,ipmP->Inf_nb,ipmP->dlamdai,&step_l,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 114 alpha = PetscMin(step_s,step_l); 115 alpha = PetscMin(alpha,1.0); 116 ipmP->alpha1 = alpha; 117 } else { 118 ipmP->alpha1 = alpha = 1.0; 119 } 120 121 122 /* x_aff = x + alpha*d */ 123 ierr = VecCopy(tao->solution,ipmP->save_x); CHKERRQ(ierr); 124 if (ipmP->me > 0) { 125 ierr = VecCopy(ipmP->lamdae,ipmP->save_lamdae); CHKERRQ(ierr); 126 } 127 if (ipmP->nb > 0) { 128 ierr = VecCopy(ipmP->lamdai,ipmP->save_lamdai); CHKERRQ(ierr); 129 ierr = VecCopy(ipmP->s,ipmP->save_s); CHKERRQ(ierr); 130 } 131 132 ierr = VecAXPY(tao->solution,alpha,tao->stepdirection); CHKERRQ(ierr); 133 if (ipmP->me > 0) { 134 ierr = VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae); CHKERRQ(ierr); 135 } 136 if (ipmP->nb > 0) { 137 ierr = VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai); CHKERRQ(ierr); 138 ierr = VecAXPY(ipmP->s,alpha,ipmP->ds); CHKERRQ(ierr); 139 } 140 141 142 /* Recompute kkt to find centering parameter sigma = (new_mu/old_mu)^3 */ 143 if (ipmP->mu == 0.0) { 144 sigma = 0.0; 145 } else { 146 sigma = 1.0/ipmP->mu; 147 } 148 ierr = IPMComputeKKT(tao); CHKERRQ(ierr); 149 sigma *= ipmP->mu; 150 sigma*=sigma*sigma; 151 152 /* revert kkt info */ 153 ierr = VecCopy(ipmP->save_x,tao->solution); CHKERRQ(ierr); 154 if (ipmP->me > 0) { 155 ierr = VecCopy(ipmP->save_lamdae,ipmP->lamdae); CHKERRQ(ierr); 156 } 157 if (ipmP->nb > 0) { 158 ierr = VecCopy(ipmP->save_lamdai,ipmP->lamdai); CHKERRQ(ierr); 159 ierr = VecCopy(ipmP->save_s,ipmP->s); CHKERRQ(ierr); 160 } 161 ierr = IPMComputeKKT(tao); CHKERRQ(ierr); 162 163 /* update rhs with new complementarity vector */ 164 if (ipmP->nb > 0) { 165 ierr = VecCopy(ipmP->complementarity,ipmP->rhs_s); CHKERRQ(ierr); 166 ierr = VecScale(ipmP->rhs_s,-1.0); CHKERRQ(ierr); 167 ierr = VecShift(ipmP->rhs_s,sigma*ipmP->mu); CHKERRQ(ierr); 168 } 169 ierr = IPMGatherRHS(tao,ipmP->bigrhs,PETSC_NULL,PETSC_NULL, 170 PETSC_NULL,ipmP->rhs_s); CHKERRQ(ierr); 171 172 /* solve K * step = rhs */ 173 ierr = KSPSetOperators(tao->ksp,ipmP->K,ipmP->K,DIFFERENT_NONZERO_PATTERN); CHKERRQ(ierr); 174 ierr = KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep);CHKERRQ(ierr); 175 176 ierr = IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds, 177 ipmP->dlamdae,ipmP->dlamdai); CHKERRQ(ierr); 178 ierr = KSPGetIterationNumber(tao->ksp,&its); CHKERRQ(ierr); 179 #if defined DEBUG_KKT2 180 PetscPrintf(PETSC_COMM_WORLD,"rhs_lamdai\n"); 181 VecView(ipmP->rhs_lamdai,0); 182 ierr = VecView(ipmP->bigrhs,0); 183 ierr = VecView(ipmP->bigstep,0); 184 #endif 185 tao->ksp_its += its; 186 187 188 if (ipmP->nb > 0) { 189 /* Get max step size and apply frac-to-boundary */ 190 tau = PetscMax(ipmP->taumin,1.0-ipmP->mu); 191 tau = PetscMin(tau,1.0); 192 if (tau != 1.0) { 193 ierr = VecScale(ipmP->s,tau); CHKERRQ(ierr); 194 ierr = VecScale(ipmP->lamdai,tau); CHKERRQ(ierr); 195 } 196 ierr = VecStepBoundInfo(ipmP->s,ipmP->Zero_nb,ipmP->Inf_nb,ipmP->ds,&step_s,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 197 ierr = VecStepBoundInfo(ipmP->lamdai,ipmP->Zero_nb,ipmP->Inf_nb,ipmP->dlamdai,&step_l,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 198 if (tau != 1.0) { 199 ierr = VecCopy(ipmP->save_s,ipmP->s); CHKERRQ(ierr); 200 ierr = VecCopy(ipmP->save_lamdai,ipmP->lamdai); CHKERRQ(ierr); 201 } 202 alpha = PetscMin(step_s,step_l); 203 alpha = PetscMin(alpha,1.0); 204 } else { 205 alpha = 1.0; 206 } 207 ipmP->alpha2 = alpha; 208 /* TODO make phi_target meaningful */ 209 phi_target = ipmP->dec * ipmP->phi; 210 for (i=0; i<11;i++) { 211 #if defined DEBUG_KKT2 212 PetscPrintf(PETSC_COMM_WORLD,"alpha2=%g\n",alpha); 213 PetscPrintf(PETSC_COMM_WORLD,"old point:\n"); 214 VecView(tao->solution,0); 215 VecView(ipmP->lamdae,0); 216 VecView(ipmP->s,0); 217 VecView(ipmP->lamdai,0); 218 #endif 219 ierr = VecAXPY(tao->solution,alpha,tao->stepdirection); CHKERRQ(ierr); 220 if (ipmP->nb > 0) { 221 ierr = VecAXPY(ipmP->s,alpha,ipmP->ds); CHKERRQ(ierr); 222 ierr = VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai); CHKERRQ(ierr); 223 } 224 if (ipmP->me > 0) { 225 ierr = VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae); CHKERRQ(ierr); 226 } 227 #if defined DEBUG_KKT 228 PetscPrintf(PETSC_COMM_WORLD,"step direction:\n"); 229 VecView(tao->stepdirection,0); 230 /* VecView(ipmP->dlamdae,0); 231 VecView(ipmP->ds,0); 232 VecView(ipmP->dlamdai,0); 233 234 PetscPrintf(PETSC_COMM_WORLD,"New iterate:\n"); 235 VecView(tao->solution,0); 236 VecView(ipmP->lamdae,0); 237 VecView(ipmP->s,0); 238 VecView(ipmP->lamdai,0); */ 239 #endif 240 /* update dual variables */ 241 if (ipmP->me > 0) { 242 ierr = VecCopy(ipmP->lamdae,tao->DE); CHKERRQ(ierr); 243 } 244 /* TODO: fix 245 if (ipmP->nb > 0) { 246 ierr = VecScatterBegin 247 PetscInt lstart,lend; 248 249 ierr = VecGetOwnershipRange(ipmP->lamdai,&lstart,&lend); 250 ierr = VecGetArray(ipmP->lamdai,&li); CHKERRQ(ierr); 251 ierr = VecGetArray(tao->DI,&di); CHKERRQ(ierr); 252 for (j=lstart;j<lend;j++) { 253 if (j < ipmP->nilb) { 254 di[j] = li[j]; 255 } 256 } 257 ierr = VecRestoreArray(ipmP->lamdai,&li); CHKERRQ(ierr); 258 ierr = VecRestoreArray(tao->DI,&di); CHKERRQ(ierr); 259 } 260 */ 261 262 263 ierr = IPMEvaluate(tao); CHKERRQ(ierr); 264 ierr = IPMComputeKKT(tao); CHKERRQ(ierr); 265 if (ipmP->phi <= phi_target) break; 266 alpha /= 2.0; 267 } 268 269 ierr = TaoMonitor(tao,iter,ipmP->kkt_f,ipmP->phi,0.0,stepsize,&reason); 270 iter++; 271 CHKERRQ(ierr); 272 273 } 274 275 PetscFunctionReturn(0); 276 } 277 278 #undef __FUNCT__ 279 #define __FUNCT__ "TaoSetup_IPM" 280 static PetscErrorCode TaoSetup_IPM(TaoSolver tao) 281 { 282 TAO_IPM *ipmP = (TAO_IPM*)tao->data; 283 PetscErrorCode ierr; 284 285 PetscFunctionBegin; 286 ipmP->nb = ipmP->mi = ipmP->me = 0; 287 ipmP->K=0; 288 ierr = VecGetSize(tao->solution,&ipmP->n); CHKERRQ(ierr); 289 if (!tao->gradient) { 290 ierr = VecDuplicate(tao->solution, &tao->gradient); CHKERRQ(ierr); 291 ierr = VecDuplicate(tao->solution, &tao->stepdirection); CHKERRQ(ierr); 292 ierr = VecDuplicate(tao->solution, &ipmP->rd); CHKERRQ(ierr); 293 ierr = VecDuplicate(tao->solution, &ipmP->rhs_x); CHKERRQ(ierr); 294 ierr = VecDuplicate(tao->solution, &ipmP->work); CHKERRQ(ierr); 295 ierr = VecDuplicate(tao->solution, &ipmP->save_x); CHKERRQ(ierr); 296 } 297 298 if (tao->constraints_equality) { 299 ierr = VecGetSize(tao->constraints_equality,&ipmP->me); CHKERRQ(ierr); 300 ierr = VecDuplicate(tao->constraints_equality,&ipmP->lamdae); CHKERRQ(ierr); 301 ierr = VecDuplicate(tao->constraints_equality,&ipmP->dlamdae); CHKERRQ(ierr); 302 ierr = VecDuplicate(tao->constraints_equality,&ipmP->rhs_lamdae); CHKERRQ(ierr); 303 ierr = VecDuplicate(tao->constraints_equality,&ipmP->save_lamdae); CHKERRQ(ierr); 304 ierr = VecDuplicate(tao->constraints_equality,&ipmP->rpe); CHKERRQ(ierr); 305 ierr = VecDuplicate(tao->constraints_equality,&tao->DE); CHKERRQ(ierr); 306 } 307 if (tao->constraints_inequality) { 308 ierr = VecDuplicate(tao->constraints_inequality,&tao->DI); CHKERRQ(ierr); 309 } 310 311 PetscFunctionReturn(0); 312 } 313 314 #undef __FUNCT__ 315 #define __FUNCT__ "IPMInitializeBounds" 316 static PetscErrorCode IPMInitializeBounds(TaoSolver tao) 317 { 318 TAO_IPM *ipmP = (TAO_IPM*)tao->data; 319 Vec xtmp; 320 /* PetscInt cstart,cend; */ /* ci (userci + lb + ub) */ 321 PetscInt xstart,xend; 322 PetscInt ucstart,ucend; /* user ci */ 323 PetscInt ucestart,uceend; /* user ce */ 324 PetscInt sstart,send; 325 PetscInt bigsize; 326 PetscInt i,counter,nloc; 327 PetscInt *cind,*xind,*ucind,*uceind,*stepind; 328 VecType vtype; 329 const PetscInt *xli,*xui; 330 PetscInt xl_offset,xu_offset; 331 IS bigxl,bigxu,isuc,isc,isx,sis,is1; 332 PetscErrorCode ierr; 333 334 MPI_Comm comm; 335 PetscFunctionBegin; 336 cind=xind=ucind=uceind=stepind=0; 337 ipmP->mi=0; 338 ipmP->nxlb=0; 339 ipmP->nxub=0; 340 ipmP->nb=0; 341 ipmP->nslack=0; 342 343 ierr = VecDuplicate(tao->solution,&xtmp); CHKERRQ(ierr); 344 if (!tao->XL && !tao->XU && tao->ops->computebounds) { 345 ierr = TaoComputeVariableBounds(tao); CHKERRQ(ierr); 346 } 347 if (tao->XL) { 348 ierr = VecSet(xtmp,TAO_NINFINITY); CHKERRQ(ierr); 349 ierr = VecWhichGreaterThan(tao->XL,xtmp,&ipmP->isxl); CHKERRQ(ierr); 350 ierr = ISGetSize(ipmP->isxl,&ipmP->nxlb); CHKERRQ(ierr); 351 } else { 352 ipmP->nxlb=0; 353 } 354 if (tao->XU) { 355 ierr = VecSet(xtmp,TAO_INFINITY); CHKERRQ(ierr); 356 ierr = VecWhichLessThan(tao->XU,xtmp,&ipmP->isxu); CHKERRQ(ierr); 357 ierr = ISGetSize(ipmP->isxu,&ipmP->nxub); CHKERRQ(ierr); 358 } else { 359 ipmP->nxub=0; 360 } 361 ierr = VecDestroy(&xtmp); CHKERRQ(ierr); 362 if (tao->constraints_inequality) { 363 ierr = VecGetSize(tao->constraints_inequality,&ipmP->mi); CHKERRQ(ierr); 364 } else { 365 ipmP->mi = 0; 366 } 367 #if defined DEBUG_K 368 PetscPrintf(PETSC_COMM_WORLD,"isxl:\n"); 369 if (ipmP->nxlb) { 370 ISView(ipmP->isxl,0); 371 } 372 PetscPrintf(PETSC_COMM_WORLD,"isxu:\n"); 373 if (ipmP->nxub) { 374 ISView(ipmP->isxu,0); 375 } 376 377 #endif 378 ipmP->nb = ipmP->nxlb + ipmP->nxub + ipmP->mi; 379 380 comm = ((PetscObject)(tao->solution))->comm; 381 382 bigsize = ipmP->n+2*ipmP->nb+ipmP->me; 383 ierr = PetscMalloc(bigsize*sizeof(PetscInt),&stepind); CHKERRQ(ierr); 384 ierr = PetscMalloc(ipmP->n*sizeof(PetscInt),&xind); CHKERRQ(ierr); 385 ierr = PetscMalloc(ipmP->me*sizeof(PetscInt),&uceind); CHKERRQ(ierr); 386 ierr = VecGetOwnershipRange(tao->solution,&xstart,&xend);CHKERRQ(ierr); 387 388 if (ipmP->nb > 0) { 389 ierr = VecCreate(comm,&ipmP->s); CHKERRQ(ierr); 390 ierr = VecSetSizes(ipmP->s,PETSC_DECIDE,ipmP->nb); CHKERRQ(ierr); 391 ierr = VecSetFromOptions(ipmP->s); CHKERRQ(ierr); 392 ierr = VecDuplicate(ipmP->s,&ipmP->ds); CHKERRQ(ierr); 393 ierr = VecDuplicate(ipmP->s,&ipmP->rhs_s); CHKERRQ(ierr); 394 ierr = VecDuplicate(ipmP->s,&ipmP->complementarity); CHKERRQ(ierr); 395 ierr = VecDuplicate(ipmP->s,&ipmP->ci); CHKERRQ(ierr); 396 397 ierr = VecDuplicate(ipmP->s,&ipmP->lamdai); CHKERRQ(ierr); 398 ierr = VecDuplicate(ipmP->s,&ipmP->dlamdai); CHKERRQ(ierr); 399 ierr = VecDuplicate(ipmP->s,&ipmP->rhs_lamdai); CHKERRQ(ierr); 400 ierr = VecDuplicate(ipmP->s,&ipmP->save_lamdai); CHKERRQ(ierr); 401 402 403 ierr = VecDuplicate(ipmP->s,&ipmP->save_s); CHKERRQ(ierr); 404 ierr = VecDuplicate(ipmP->s,&ipmP->rpi); CHKERRQ(ierr); 405 ierr = VecDuplicate(ipmP->s,&ipmP->Zero_nb); CHKERRQ(ierr); 406 ierr = VecSet(ipmP->Zero_nb,0.0); CHKERRQ(ierr); 407 ierr = VecDuplicate(ipmP->s,&ipmP->One_nb); CHKERRQ(ierr); 408 ierr = VecSet(ipmP->One_nb,1.0); CHKERRQ(ierr); 409 ierr = VecDuplicate(ipmP->s,&ipmP->Inf_nb); CHKERRQ(ierr); 410 ierr = VecSet(ipmP->Inf_nb,TAO_INFINITY); CHKERRQ(ierr); 411 412 ierr = PetscMalloc(ipmP->nb*sizeof(PetscInt),&cind); CHKERRQ(ierr); 413 ierr = PetscMalloc(ipmP->mi*sizeof(PetscInt),&ucind); CHKERRQ(ierr); 414 ierr = VecGetOwnershipRange(ipmP->s,&sstart,&send);CHKERRQ(ierr); 415 416 417 418 if (ipmP->mi > 0) { 419 ierr = VecGetOwnershipRange(tao->constraints_inequality,&ucstart,&ucend); CHKERRQ(ierr); 420 counter=0; 421 for (i=ucstart;i<ucend;i++) { 422 cind[counter++] = i; 423 } 424 ierr = ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isuc); CHKERRQ(ierr); 425 ierr = ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isc); CHKERRQ(ierr); 426 ierr = VecScatterCreate(tao->constraints_inequality,isuc,ipmP->ci,isc,&ipmP->ci_scat); CHKERRQ(ierr); 427 428 ierr = ISDestroy(&isuc); 429 ierr = ISDestroy(&isc); 430 } 431 /* need to know how may xbound indices are on each process */ 432 /* TODO better way */ 433 if (ipmP->nxlb) { 434 ierr = ISAllGather(ipmP->isxl,&bigxl);CHKERRQ(ierr); 435 ierr = ISGetIndices(bigxl,&xli);CHKERRQ(ierr); 436 /* find offsets for this processor */ 437 xl_offset = ipmP->mi; 438 for (i=0;i<ipmP->nxlb;i++) { 439 if (xli[i] < xstart) { 440 xl_offset++; 441 } else break; 442 } 443 ierr = ISRestoreIndices(bigxl,&xli);CHKERRQ(ierr); 444 445 ierr = ISGetIndices(ipmP->isxl,&xli);CHKERRQ(ierr); 446 ierr = ISGetLocalSize(ipmP->isxl,&nloc);CHKERRQ(ierr); 447 for (i=0;i<nloc;i++) { 448 xind[i] = xli[i]; 449 cind[i] = xl_offset+i; 450 } 451 452 ierr = ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx); CHKERRQ(ierr); 453 ierr = ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr); 454 ierr = VecScatterCreate(tao->XL,isx,ipmP->ci,isc,&ipmP->xl_scat); CHKERRQ(ierr); 455 ierr = ISDestroy(&isx);CHKERRQ(ierr); 456 ierr = ISDestroy(&isc);CHKERRQ(ierr); 457 ierr = ISDestroy(&bigxl);CHKERRQ(ierr); 458 } 459 460 461 if (ipmP->nxub) { 462 ierr = ISAllGather(ipmP->isxu,&bigxu);CHKERRQ(ierr); 463 ierr = ISGetIndices(bigxu,&xui);CHKERRQ(ierr); 464 /* find offsets for this processor */ 465 xu_offset = ipmP->mi + ipmP->nxlb; 466 for (i=0;i<ipmP->nxub;i++) { 467 if (xui[i] < xstart) { 468 xu_offset++; 469 } else break; 470 } 471 ierr = ISRestoreIndices(bigxu,&xui);CHKERRQ(ierr); 472 473 ierr = ISGetIndices(ipmP->isxu,&xui);CHKERRQ(ierr); 474 ierr = ISGetLocalSize(ipmP->isxu,&nloc);CHKERRQ(ierr); 475 for (i=0;i<nloc;i++) { 476 xind[i] = xui[i]; 477 cind[i] = xu_offset+i; 478 } 479 480 ierr = ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx); CHKERRQ(ierr); 481 ierr = ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr); 482 ierr = VecScatterCreate(tao->XU,isx,ipmP->ci,isc,&ipmP->xu_scat); CHKERRQ(ierr); 483 ierr = ISDestroy(&isx);CHKERRQ(ierr); 484 ierr = ISDestroy(&isc);CHKERRQ(ierr); 485 ierr = ISDestroy(&bigxu);CHKERRQ(ierr); 486 } 487 } 488 489 490 ierr = VecCreate(comm,&ipmP->bigrhs); CHKERRQ(ierr); 491 ierr = VecGetType(tao->solution,&vtype); CHKERRQ(ierr); 492 ierr = VecSetType(ipmP->bigrhs,vtype); CHKERRQ(ierr); 493 ierr = VecSetSizes(ipmP->bigrhs,PETSC_DECIDE,bigsize); CHKERRQ(ierr); 494 ierr = VecSetFromOptions(ipmP->bigrhs); CHKERRQ(ierr); 495 ierr = VecDuplicate(ipmP->bigrhs,&ipmP->bigstep); CHKERRQ(ierr); 496 497 /* create scatters for step->x and x->rhs */ 498 for (i=xstart;i<xend;i++) { 499 stepind[i-xstart] = i; 500 xind[i-xstart] = i; 501 } 502 ierr = ISCreateGeneral(comm,xend-xstart,stepind,PETSC_COPY_VALUES,&sis); CHKERRQ(ierr); 503 ierr = ISCreateGeneral(comm,xend-xstart,xind,PETSC_COPY_VALUES,&is1); CHKERRQ(ierr); 504 ierr = VecScatterCreate(ipmP->bigstep,sis,tao->solution,is1,&ipmP->step1); CHKERRQ(ierr); 505 ierr = VecScatterCreate(tao->solution,is1,ipmP->bigrhs,sis,&ipmP->rhs1); CHKERRQ(ierr); 506 ierr = ISDestroy(&sis); CHKERRQ(ierr); 507 ierr = ISDestroy(&is1); CHKERRQ(ierr); 508 509 if (ipmP->nb > 0) { 510 for (i=sstart;i<send;i++) { 511 stepind[i-sstart] = i+ipmP->n; 512 cind[i-sstart] = i; 513 } 514 ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis); CHKERRQ(ierr); 515 ierr = ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1); CHKERRQ(ierr); 516 ierr = VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step2); CHKERRQ(ierr); 517 ierr = ISDestroy(&sis); CHKERRQ(ierr); 518 519 for (i=sstart;i<send;i++) { 520 stepind[i-sstart] = i+ipmP->n+ipmP->me; 521 cind[i-sstart] = i; 522 } 523 ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis); CHKERRQ(ierr); 524 ierr = VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs3); CHKERRQ(ierr); 525 ierr = ISDestroy(&sis); CHKERRQ(ierr); 526 ierr = ISDestroy(&is1); CHKERRQ(ierr); 527 528 529 } 530 531 if (ipmP->me > 0) { 532 ierr = VecGetOwnershipRange(tao->constraints_equality,&ucestart,&uceend); CHKERRQ(ierr); 533 for (i=ucestart;i<uceend;i++) { 534 stepind[i-ucestart] = i + ipmP->n+ipmP->nb; 535 uceind[i-ucestart] = i; 536 } 537 538 ierr = ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis); CHKERRQ(ierr); 539 ierr = ISCreateGeneral(comm,uceend-ucestart,uceind,PETSC_COPY_VALUES,&is1); CHKERRQ(ierr); 540 ierr = VecScatterCreate(ipmP->bigstep,sis,tao->constraints_equality,is1,&ipmP->step3); CHKERRQ(ierr); 541 ierr = ISDestroy(&sis); CHKERRQ(ierr); 542 543 544 for (i=ucestart;i<uceend;i++) { 545 stepind[i-ucestart] = i + ipmP->n; 546 } 547 548 ierr = ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis); CHKERRQ(ierr); 549 ierr = VecScatterCreate(tao->constraints_equality,is1,ipmP->bigrhs,sis,&ipmP->rhs2); CHKERRQ(ierr); 550 ierr = ISDestroy(&sis); CHKERRQ(ierr); 551 ierr = ISDestroy(&is1); CHKERRQ(ierr); 552 } 553 554 if (ipmP->nb > 0) { 555 for (i=sstart;i<send;i++) { 556 stepind[i-sstart] = i + ipmP->n + ipmP->nb + ipmP->me; 557 cind[i-sstart] = i; 558 } 559 ierr = ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1); 560 ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis); CHKERRQ(ierr); 561 ierr = VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step4); CHKERRQ(ierr); 562 ierr = VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs4); CHKERRQ(ierr); 563 ierr = ISDestroy(&sis); CHKERRQ(ierr); 564 ierr = ISDestroy(&is1); CHKERRQ(ierr); 565 } 566 567 ierr = PetscFree(stepind); CHKERRQ(ierr); 568 ierr = PetscFree(cind); CHKERRQ(ierr); 569 ierr = PetscFree(ucind); CHKERRQ(ierr); 570 ierr = PetscFree(uceind); CHKERRQ(ierr); 571 ierr = PetscFree(xind); CHKERRQ(ierr); 572 573 PetscFunctionReturn(0); 574 } 575 576 #undef __FUNCT__ 577 #define __FUNCT__ "TaoDestroy_IPM" 578 static PetscErrorCode TaoDestroy_IPM(TaoSolver tao) 579 { 580 TAO_IPM *ipmP = (TAO_IPM*)tao->data; 581 PetscErrorCode ierr; 582 PetscFunctionBegin; 583 ierr = VecDestroy(&ipmP->rd); CHKERRQ(ierr); 584 ierr = VecDestroy(&ipmP->rpe); CHKERRQ(ierr); 585 ierr = VecDestroy(&ipmP->rpi); CHKERRQ(ierr); 586 ierr = VecDestroy(&ipmP->work); CHKERRQ(ierr); 587 ierr = VecDestroy(&ipmP->lamdae); CHKERRQ(ierr); 588 ierr = VecDestroy(&ipmP->lamdai); CHKERRQ(ierr); 589 ierr = VecDestroy(&ipmP->s); CHKERRQ(ierr); 590 ierr = VecDestroy(&ipmP->ds); CHKERRQ(ierr); 591 ierr = VecDestroy(&ipmP->ci); CHKERRQ(ierr); 592 593 ierr = VecDestroy(&ipmP->rhs_x); CHKERRQ(ierr); 594 ierr = VecDestroy(&ipmP->rhs_lamdae); CHKERRQ(ierr); 595 ierr = VecDestroy(&ipmP->rhs_lamdai); CHKERRQ(ierr); 596 ierr = VecDestroy(&ipmP->rhs_s); CHKERRQ(ierr); 597 598 ierr = VecDestroy(&ipmP->save_x); CHKERRQ(ierr); 599 ierr = VecDestroy(&ipmP->save_lamdae); CHKERRQ(ierr); 600 ierr = VecDestroy(&ipmP->save_lamdai); CHKERRQ(ierr); 601 ierr = VecDestroy(&ipmP->save_s); CHKERRQ(ierr); 602 603 ierr = VecScatterDestroy(&ipmP->step1); CHKERRQ(ierr); 604 ierr = VecScatterDestroy(&ipmP->step2); CHKERRQ(ierr); 605 ierr = VecScatterDestroy(&ipmP->step3); CHKERRQ(ierr); 606 ierr = VecScatterDestroy(&ipmP->step4); CHKERRQ(ierr); 607 608 ierr = VecScatterDestroy(&ipmP->rhs1); CHKERRQ(ierr); 609 ierr = VecScatterDestroy(&ipmP->rhs2); CHKERRQ(ierr); 610 ierr = VecScatterDestroy(&ipmP->rhs3); CHKERRQ(ierr); 611 ierr = VecScatterDestroy(&ipmP->rhs4); CHKERRQ(ierr); 612 613 ierr = VecScatterDestroy(&ipmP->ci_scat); CHKERRQ(ierr); 614 ierr = VecScatterDestroy(&ipmP->xl_scat); CHKERRQ(ierr); 615 ierr = VecScatterDestroy(&ipmP->xu_scat); CHKERRQ(ierr); 616 617 ierr = VecDestroy(&ipmP->dlamdai); CHKERRQ(ierr); 618 ierr = VecDestroy(&ipmP->dlamdae); CHKERRQ(ierr); 619 ierr = VecDestroy(&ipmP->Zero_nb); CHKERRQ(ierr); 620 ierr = VecDestroy(&ipmP->One_nb); CHKERRQ(ierr); 621 ierr = VecDestroy(&ipmP->Inf_nb); CHKERRQ(ierr); 622 ierr = VecDestroy(&ipmP->complementarity); CHKERRQ(ierr); 623 624 625 ierr = VecDestroy(&ipmP->bigrhs); CHKERRQ(ierr); 626 ierr = VecDestroy(&ipmP->bigstep); CHKERRQ(ierr); 627 ierr = MatDestroy(&ipmP->Ai); CHKERRQ(ierr); 628 ierr = MatDestroy(&ipmP->K); CHKERRQ(ierr); 629 ierr = ISDestroy(&ipmP->isxu); CHKERRQ(ierr); 630 ierr = ISDestroy(&ipmP->isxl); CHKERRQ(ierr); 631 ierr = PetscFree(tao->data); CHKERRQ(ierr); 632 tao->data = PETSC_NULL; 633 PetscFunctionReturn(0); 634 } 635 636 #undef __FUNCT__ 637 #define __FUNCT__ "TaoSetFromOptions_IPM" 638 static PetscErrorCode TaoSetFromOptions_IPM(TaoSolver tao) 639 { 640 TAO_IPM *ipmP = (TAO_IPM*)tao->data; 641 PetscErrorCode ierr; 642 PetscBool flg; 643 PetscFunctionBegin; 644 ierr = PetscOptionsHead("IPM method for constrained optimization"); CHKERRQ(ierr); 645 ierr = PetscOptionsBool("-ipm_monitorkkt","monitor kkt status",PETSC_NULL,ipmP->monitorkkt,&ipmP->monitorkkt,&flg); CHKERRQ(ierr); 646 ierr = PetscOptionsReal("-ipm_pushs","parameter to push initial slack variables away from bounds",PETSC_NULL,ipmP->pushs,&ipmP->pushs,&flg); 647 ierr = PetscOptionsReal("-ipm_pushnu","parameter to push initial (inequality) dual variables away from bounds",PETSC_NULL,ipmP->pushnu,&ipmP->pushnu,&flg); 648 ierr = PetscOptionsTail(); CHKERRQ(ierr); 649 ierr =KSPSetFromOptions(tao->ksp); CHKERRQ(ierr); 650 PetscFunctionReturn(0); 651 } 652 653 #undef __FUNCT__ 654 #define __FUNCT__ "TaoView_IPM" 655 static PetscErrorCode TaoView_IPM(TaoSolver tao, PetscViewer viewer) 656 { 657 return 0; 658 } 659 660 /* IPMObjectiveAndGradient() 661 f = d'x + 0.5 * x' * H * x 662 rd = H*x + d + Ae'*lame - Ai'*lami 663 rpe = Ae*x - be 664 rpi = Ai*x - yi - bi 665 mu = yi' * lami/mi; 666 com = yi.*lami 667 668 phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 669 */ 670 /* 671 #undef __FUNCT__ 672 #define __FUNCT__ "IPMObjective" 673 static PetscErrorCode IPMObjective(TaoLineSearch ls, Vec X, PetscReal *f, void *tptr) 674 { 675 TaoSolver tao = (TaoSolver)tptr; 676 TAO_IPM *ipmP = (TAO_IPM*)tao->data; 677 PetscErrorCode ierr; 678 PetscFunctionBegin; 679 ierr = IPMComputeKKT(tao); CHKERRQ(ierr); 680 *f = ipmP->phi; 681 PetscFunctionReturn(0); 682 } 683 */ 684 685 /* 686 f = d'x + 0.5 * x' * H * x 687 rd = H*x + d + Ae'*lame - Ai'*lami 688 Ai = jac_ineq 689 I (w/lb) 690 -I (w/ub) 691 692 rpe = ce 693 rpi = ci - s; 694 com = s.*lami 695 mu = yi' * lami/mi; 696 697 phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 698 */ 699 #undef __FUNCT__ 700 #define __FUNCT__ "IPMComputeKKT" 701 static PetscErrorCode IPMComputeKKT(TaoSolver tao) 702 { 703 TAO_IPM *ipmP = (TAO_IPM *)tao->data; 704 PetscScalar norm; 705 PetscErrorCode ierr; 706 ierr = VecCopy(tao->gradient,ipmP->rd); CHKERRQ(ierr); 707 708 if (ipmP->me > 0) { 709 /* rd = gradient + Ae'*lamdae */ 710 ierr = MatMultTranspose(tao->jacobian_equality,ipmP->lamdae,ipmP->work); CHKERRQ(ierr); 711 ierr = VecAXPY(ipmP->rd, 1.0, ipmP->work); CHKERRQ(ierr); 712 713 #if defined DEBUG_KKT 714 PetscPrintf(PETSC_COMM_WORLD,"\nAe.lamdae\n"); 715 ierr = VecView(ipmP->work,0); 716 #endif 717 /* rpe = ce(x) */ 718 ierr = VecCopy(tao->constraints_equality,ipmP->rpe); CHKERRQ(ierr); 719 720 } 721 if (ipmP->nb > 0) { 722 /* rd = rd - Ai'*lamdai */ 723 ierr = MatMultTranspose(ipmP->Ai,ipmP->lamdai,ipmP->work); CHKERRQ(ierr); 724 ierr = VecAXPY(ipmP->rd, -1.0, ipmP->work); CHKERRQ(ierr); 725 #if defined DEBUG_KKT 726 PetscPrintf(PETSC_COMM_WORLD,"\nAi\n"); 727 ierr = MatView(ipmP->Ai,0); 728 PetscPrintf(PETSC_COMM_WORLD,"\nAi.lamdai\n"); 729 ierr = VecView(ipmP->work,0); 730 #endif 731 /* rpi = cin - s */ 732 ierr = VecCopy(ipmP->ci,ipmP->rpi); CHKERRQ(ierr); 733 ierr = VecAXPY(ipmP->rpi, -1.0, ipmP->s); CHKERRQ(ierr); 734 735 /* com = s .* lami */ 736 ierr = VecPointwiseMult(ipmP->complementarity, ipmP->s,ipmP->lamdai); CHKERRQ(ierr); 737 738 } 739 /* phi = ||rd; rpe; rpi; com|| */ 740 ierr = VecDot(ipmP->rd,ipmP->rd,&norm); CHKERRQ(ierr); 741 ipmP->phi = norm; 742 if (ipmP->me > 0 ) { 743 ierr = VecDot(ipmP->rpe,ipmP->rpe,&norm); CHKERRQ(ierr); 744 ipmP->phi += norm; 745 } 746 if (ipmP->nb > 0) { 747 ierr = VecDot(ipmP->rpi,ipmP->rpi,&norm); CHKERRQ(ierr); 748 ipmP->phi += norm; 749 ierr = VecDot(ipmP->complementarity,ipmP->complementarity,&norm); CHKERRQ(ierr); 750 ipmP->phi += norm; 751 /* mu = s'*lami/nb */ 752 ierr = VecDot(ipmP->s,ipmP->lamdai,&ipmP->mu); CHKERRQ(ierr); 753 ipmP->mu /= ipmP->nb; 754 } else { 755 ipmP->mu = 1.0; 756 } 757 758 ipmP->phi = PetscSqrtScalar(ipmP->phi); 759 #if defined DEBUG_KKT 760 if (ipmP->monitorkkt) { 761 ierr = PetscPrintf(PETSC_COMM_WORLD,"obj=%G,\tphi = %G,\tmu=%G\talpha1=%G\talpha2=%G\n",ipmP->kkt_f,ipmP->phi,ipmP->mu,ipmP->alpha1,ipmP->alpha2); 762 } 763 CHKMEMQ; 764 PetscPrintf(PETSC_COMM_WORLD,"\ngradient\n"); 765 ierr = VecView(tao->gradient,0); 766 PetscPrintf(PETSC_COMM_WORLD,"\nrd\n"); 767 ierr = VecView(ipmP->rd,0); 768 PetscPrintf(PETSC_COMM_WORLD,"\nrpe\n"); 769 ierr = VecView(ipmP->rpe,0); 770 PetscPrintf(PETSC_COMM_WORLD,"\nrpi\n"); 771 ierr = VecView(ipmP->rpi,0); 772 PetscPrintf(PETSC_COMM_WORLD,"\ncomplementarity\n"); 773 ierr = VecView(ipmP->complementarity,0); 774 #endif 775 PetscFunctionReturn(0); 776 } 777 778 #undef __FUNCT__ 779 #define __FUNCT__ "IPMEvaluate" 780 /* evaluate user info at current point */ 781 PetscErrorCode IPMEvaluate(TaoSolver tao) 782 { 783 TAO_IPM *ipmP = (TAO_IPM *)tao->data; 784 PetscErrorCode ierr; 785 PetscFunctionBegin; 786 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&ipmP->kkt_f,tao->gradient); CHKERRQ(ierr); 787 ierr = TaoComputeHessian(tao,tao->solution,&tao->hessian,&tao->hessian_pre,&ipmP->Hflag); CHKERRQ(ierr); 788 789 if (ipmP->me > 0) { 790 ierr = TaoComputeEqualityConstraints(tao,tao->solution,tao->constraints_equality); 791 ierr = TaoComputeJacobianEquality(tao,tao->solution,&tao->jacobian_equality,&tao->jacobian_equality_pre,&ipmP->Aiflag); CHKERRQ(ierr); 792 } 793 if (ipmP->mi > 0) { 794 ierr = TaoComputeInequalityConstraints(tao,tao->solution,tao->constraints_inequality); 795 ierr = TaoComputeJacobianInequality(tao,tao->solution,&tao->jacobian_inequality,&tao->jacobian_inequality_pre,&ipmP->Aeflag); CHKERRQ(ierr); 796 797 } 798 if (ipmP->nb > 0) { 799 /* Ai' = jac_ineq | I (w/lb) | -I (w/ub) */ 800 ierr = IPMUpdateAi(tao); CHKERRQ(ierr); 801 } 802 PetscFunctionReturn(0); 803 } 804 805 #undef __FUNCT__ 806 #define __FUNCT__ "IPMPushInitialPoint" 807 /* Push initial point away from bounds */ 808 PetscErrorCode IPMPushInitialPoint(TaoSolver tao) 809 { 810 TAO_IPM *ipmP = (TAO_IPM *)tao->data; 811 PetscErrorCode ierr; 812 PetscFunctionBegin; 813 814 ierr = TaoComputeVariableBounds(tao); CHKERRQ(ierr); 815 if (tao->XL && tao->XU) { 816 ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution); CHKERRQ(ierr); 817 } 818 if (ipmP->nb > 0) { 819 ierr = VecSet(ipmP->s,ipmP->pushs); CHKERRQ(ierr); 820 ierr = VecSet(ipmP->lamdai,ipmP->pushnu); CHKERRQ(ierr); 821 if (ipmP->mi > 0) { 822 ierr = VecSet(tao->DI,ipmP->pushnu); CHKERRQ(ierr); 823 } 824 } 825 if (ipmP->me > 0) { 826 ierr = VecSet(tao->DE,1.0); CHKERRQ(ierr); 827 ierr = VecSet(ipmP->lamdae,1.0); CHKERRQ(ierr); 828 } 829 830 831 PetscFunctionReturn(0); 832 } 833 834 #undef __FUNCT__ 835 #define __FUNCT__ "IPMUpdateAi" 836 PetscErrorCode IPMUpdateAi(TaoSolver tao) 837 { 838 /* Ai = Ji 839 I (w/lb) 840 -I (w/ub) */ 841 842 /* Ci = user->ci 843 Xi - lb (w/lb) 844 -Xi + ub (w/ub) */ 845 846 TAO_IPM *ipmP = (TAO_IPM *)tao->data; 847 MPI_Comm comm; 848 PetscInt i; 849 PetscScalar newval; 850 PetscInt newrow,newcol,ncols; 851 const PetscScalar *vals; 852 const PetscInt *cols; 853 PetscInt astart,aend,jstart,jend; 854 PetscInt *nonzeros; 855 PetscInt r2,r3,r4; 856 PetscMPIInt mpisize; 857 PetscErrorCode ierr; 858 859 PetscFunctionBegin; 860 CHKMEMQ; 861 r2 = ipmP->mi; 862 r3 = r2 + ipmP->nxlb; 863 r4 = r3 + ipmP->nxub; 864 865 if (!ipmP->nb) { 866 PetscFunctionReturn(0); 867 } 868 CHKMEMQ; 869 870 /* Create Ai matrix if it doesn't exist yet */ 871 if (!ipmP->Ai) { 872 comm = ((PetscObject)(tao->solution))->comm; 873 ierr = PetscMalloc(ipmP->nb*sizeof(PetscInt),&nonzeros); CHKERRQ(ierr); 874 ierr = MPI_Comm_size(comm,&mpisize); 875 if (mpisize == 1) { 876 for (i=0;i<ipmP->mi;i++) { 877 ierr = MatGetRow(tao->jacobian_inequality,i,&ncols,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 878 nonzeros[i] = ncols; 879 ierr = MatRestoreRow(tao->jacobian_inequality,i,&ncols,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 880 } 881 for (i=r2;i<r4;i++) { 882 nonzeros[i] = 1; 883 } 884 } 885 ierr = MatCreate(comm,&ipmP->Ai); CHKERRQ(ierr); 886 ierr = MatSetType(ipmP->Ai,MATAIJ); CHKERRQ(ierr); 887 ierr = MatSetSizes(ipmP->Ai,PETSC_DECIDE,PETSC_DECIDE,ipmP->nb,ipmP->n);CHKERRQ(ierr); 888 ierr = MatSetFromOptions(ipmP->Ai); CHKERRQ(ierr); 889 ierr = MatMPIAIJSetPreallocation(ipmP->Ai,ipmP->nb,PETSC_NULL,ipmP->nb,PETSC_NULL); 890 ierr = MatSeqAIJSetPreallocation(ipmP->Ai,PETSC_DEFAULT,nonzeros);CHKERRQ(ierr); 891 if (mpisize ==1) { 892 ierr = PetscFree(nonzeros);CHKERRQ(ierr); 893 } 894 } 895 896 897 /* Copy values from user jacobian to Ai */ 898 ierr = MatGetOwnershipRange(ipmP->Ai,&astart,&aend); CHKERRQ(ierr); 899 900 /* Ai w/lb */ 901 if (ipmP->mi) { 902 ierr = MatZeroEntries(ipmP->Ai); CHKERRQ(ierr); 903 ierr = MatGetOwnershipRange(tao->jacobian_inequality,&jstart,&jend); CHKERRQ(ierr); 904 for (i=jstart;i<jend;i++) { 905 ierr = MatGetRow(tao->jacobian_inequality,i,&ncols,&cols,&vals); CHKERRQ(ierr); 906 newrow = i; 907 ierr = MatSetValues(ipmP->Ai,1,&newrow,ncols,cols,vals,INSERT_VALUES); CHKERRQ(ierr); 908 ierr = MatRestoreRow(tao->jacobian_inequality,i,&ncols,&cols,&vals); CHKERRQ(ierr); 909 } 910 } 911 912 913 /* I w/ xlb */ 914 if (ipmP->nxlb) { 915 for (i=0;i<ipmP->nxlb;i++) { 916 if (i>=astart && i<aend) { 917 newrow = i+r2; 918 newcol = i; 919 newval = 1.0; 920 ierr = MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES); CHKERRQ(ierr); 921 } 922 } 923 } 924 if (ipmP->nxub) { 925 /* I w/ xub */ 926 for (i=0;i<ipmP->nxub;i++) { 927 if (i>=astart && i<aend) { 928 newrow = i+r3; 929 newcol = i; 930 newval = -1.0; 931 ierr = MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES); CHKERRQ(ierr); 932 } 933 } 934 } 935 936 937 ierr = MatAssemblyBegin(ipmP->Ai,MAT_FINAL_ASSEMBLY); 938 ierr = MatAssemblyEnd(ipmP->Ai,MAT_FINAL_ASSEMBLY); 939 CHKMEMQ; 940 941 ierr = VecSet(ipmP->ci,0.0); CHKERRQ(ierr); 942 943 /* user ci */ 944 if (ipmP->mi > 0) { 945 ierr = VecScatterBegin(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr); 946 ierr = VecScatterEnd(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr); 947 } 948 if (!ipmP->work){ 949 VecDuplicate(tao->solution,&ipmP->work); 950 } 951 ierr = VecCopy(tao->solution,ipmP->work);CHKERRQ(ierr); 952 if (tao->XL) { 953 ierr = VecAXPY(ipmP->work,-1.0,tao->XL);CHKERRQ(ierr); 954 955 /* lower bounds on variables */ 956 if (ipmP->nxlb > 0) { 957 ierr = VecScatterBegin(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr); 958 ierr = VecScatterEnd(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr); 959 } 960 } 961 if (tao->XU) { 962 /* upper bounds on variables */ 963 ierr = VecCopy(tao->solution,ipmP->work);CHKERRQ(ierr); 964 ierr = VecScale(ipmP->work,-1.0);CHKERRQ(ierr); 965 ierr = VecAXPY(ipmP->work,1.0,tao->XU);CHKERRQ(ierr); 966 if (ipmP->nxub > 0) { 967 ierr = VecScatterBegin(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr); 968 ierr = VecScatterEnd(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr); 969 } 970 } 971 972 PetscFunctionReturn(0); 973 } 974 975 #undef __FUNCT__ 976 #define __FUNCT__ "IPMUpdateK" 977 /* create K = [ Hlag , 0 , Ae', -Ai']; 978 [Ae , 0, 0 , 0]; 979 [Ai ,-I, 0 , 0]; 980 [ 0 , S , 0, Y ]; */ 981 PetscErrorCode IPMUpdateK(TaoSolver tao) 982 { 983 TAO_IPM *ipmP = (TAO_IPM *)tao->data; 984 MPI_Comm comm; 985 PetscMPIInt mpisize; 986 PetscErrorCode ierr; 987 PetscInt i,j,row; 988 PetscInt ncols,newcol,newcols[2],newrow; 989 const PetscInt *cols; 990 const PetscReal *vals; 991 PetscReal *l,*y; 992 PetscReal *newvals; 993 PetscReal newval; 994 PetscInt subsize; 995 const PetscInt *indices; 996 PetscInt *nonzeros,*d_nonzeros,*o_nonzeros; 997 PetscInt bigsize; 998 PetscInt r1,r2,r3; 999 PetscInt c1,c2,c3; 1000 PetscInt klocalsize; 1001 PetscInt hstart,hend,kstart,kend; 1002 PetscInt aistart,aiend,aestart,aeend; 1003 PetscInt sstart,send; 1004 PetscFunctionBegin; 1005 1006 comm = ((PetscObject)(tao->solution))->comm; 1007 ierr = MPI_Comm_size(comm,&mpisize); 1008 ierr = IPMUpdateAi(tao); CHKERRQ(ierr); 1009 #if defined DEBUG_K 1010 PetscPrintf(PETSC_COMM_WORLD,"H\n"); MatView(tao->hessian,0); 1011 if (ipmP->nb) { 1012 PetscPrintf(PETSC_COMM_WORLD,"Ai\n"); MatView(ipmP->Ai,0); 1013 } 1014 if (ipmP->me) { 1015 PetscPrintf(PETSC_COMM_WORLD,"Ae\n"); MatView(tao->jacobian_equality,0); 1016 } 1017 1018 #endif 1019 /* allocate workspace */ 1020 subsize = PetscMax(ipmP->n,ipmP->nb); 1021 subsize = PetscMax(ipmP->me,subsize); 1022 subsize = PetscMax(2,subsize); 1023 ierr = PetscMalloc(sizeof(PetscInt)*subsize,&indices); CHKERRQ(ierr); 1024 ierr = PetscMalloc(sizeof(PetscReal)*subsize,&newvals); CHKERRQ(ierr); 1025 1026 r1 = c1 = ipmP->n; 1027 r2 = r1 + ipmP->me; c2 = c1 + ipmP->nb; 1028 r3 = c3 = r2 + ipmP->nb; 1029 1030 1031 1032 bigsize = ipmP->n+2*ipmP->nb+ipmP->me; 1033 ierr = VecGetOwnershipRange(ipmP->bigrhs,&kstart,&kend); CHKERRQ(ierr); 1034 ierr = MatGetOwnershipRange(tao->hessian,&hstart,&hend); CHKERRQ(ierr); 1035 klocalsize = kend-kstart; 1036 if (!ipmP->K) { 1037 if (mpisize == 1) { 1038 ierr = PetscMalloc((kend-kstart)*sizeof(PetscInt),&nonzeros);CHKERRQ(ierr); 1039 for (i=0;i<bigsize;i++) { 1040 if (i<r1) { 1041 ierr = MatGetRow(tao->hessian,i,&ncols,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 1042 nonzeros[i] = ncols; 1043 ierr = MatRestoreRow(tao->hessian,i,&ncols,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 1044 nonzeros[i] += ipmP->me+ipmP->nb; 1045 } else if (i<r2) { 1046 nonzeros[i-kstart] = ipmP->n; 1047 } else if (i<r3) { 1048 nonzeros[i-kstart] = ipmP->n+1; 1049 } else if (i<bigsize) { 1050 nonzeros[i-kstart] = 2; 1051 } 1052 } 1053 ierr = MatCreate(comm,&ipmP->K); CHKERRQ(ierr); 1054 ierr = MatSetType(ipmP->K,MATSEQAIJ); CHKERRQ(ierr); 1055 ierr = MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1056 ierr = MatSeqAIJSetPreallocation(ipmP->K,0,nonzeros); CHKERRQ(ierr); 1057 ierr = MatSetFromOptions(ipmP->K); CHKERRQ(ierr); 1058 ierr = PetscFree(nonzeros); CHKERRQ(ierr); 1059 } else { 1060 ierr = PetscMalloc((kend-kstart)*sizeof(PetscInt),&d_nonzeros); CHKERRQ(ierr); 1061 ierr = PetscMalloc((kend-kstart)*sizeof(PetscInt),&o_nonzeros); CHKERRQ(ierr); 1062 for (i=kstart;i<kend;i++) { 1063 if (i<r1) { 1064 /* TODO fix preallocation for mpi mats */ 1065 d_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,kend-kstart); 1066 o_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,bigsize-(kend-kstart)); 1067 } else if (i<r2) { 1068 d_nonzeros[i-kstart] = PetscMin(ipmP->n,kend-kstart); 1069 o_nonzeros[i-kstart] = PetscMin(ipmP->n,bigsize-(kend-kstart)); 1070 } else if (i<r3) { 1071 d_nonzeros[i-kstart] = PetscMin(ipmP->n+2,kend-kstart); 1072 o_nonzeros[i-kstart] = PetscMin(ipmP->n+2,bigsize-(kend-kstart)); 1073 } else { 1074 d_nonzeros[i-kstart] = PetscMin(2,kend-kstart); 1075 o_nonzeros[i-kstart] = PetscMin(2,bigsize-(kend-kstart)); 1076 } 1077 } 1078 ierr = MatCreate(comm,&ipmP->K); CHKERRQ(ierr); 1079 ierr = MatSetType(ipmP->K,MATMPIAIJ); CHKERRQ(ierr); 1080 ierr = MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1081 ierr = MatMPIAIJSetPreallocation(ipmP->K,0,d_nonzeros,0,o_nonzeros); CHKERRQ(ierr); 1082 ierr = PetscFree(d_nonzeros); CHKERRQ(ierr); 1083 ierr = PetscFree(o_nonzeros); CHKERRQ(ierr); 1084 ierr = MatSetFromOptions(ipmP->K); CHKERRQ(ierr); 1085 1086 } 1087 } 1088 1089 1090 ierr = MatZeroEntries(ipmP->K); CHKERRQ(ierr); 1091 /* Copy H */ 1092 for (i=hstart;i<hend;i++) { 1093 ierr = MatGetRow(tao->hessian,i,&ncols,&cols,&vals); CHKERRQ(ierr); 1094 if (ncols > 0) { 1095 ierr = MatSetValues(ipmP->K,1,&i,ncols,cols,vals,INSERT_VALUES); CHKERRQ(ierr); 1096 } 1097 ierr = MatRestoreRow(tao->hessian,i,&ncols,&cols,&vals); CHKERRQ(ierr); 1098 } 1099 1100 /* Copy Ae and Ae' */ 1101 if (ipmP->me > 0) { 1102 ierr = MatGetOwnershipRange(tao->jacobian_equality,&aestart,&aeend); CHKERRQ(ierr); 1103 for (i=aestart;i<aeend;i++) { 1104 ierr = MatGetRow(tao->jacobian_equality,i,&ncols,&cols,&vals); CHKERRQ(ierr); 1105 if (ncols > 0) { 1106 /*Ae*/ 1107 row = i+r1; 1108 ierr = MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES); CHKERRQ(ierr); 1109 /*Ae'*/ 1110 for (j=0;j<ncols;j++) { 1111 newcol = i + c2; 1112 newrow = cols[j]; 1113 newval = vals[j]; 1114 ierr = MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES); CHKERRQ(ierr); 1115 } 1116 } 1117 ierr = MatRestoreRow(tao->jacobian_equality,i,&ncols,&cols,&vals); CHKERRQ(ierr); 1118 } 1119 } 1120 1121 if (ipmP->nb > 0) { 1122 ierr = MatGetOwnershipRange(ipmP->Ai,&aistart,&aiend); CHKERRQ(ierr); 1123 /* Copy Ai,and Ai' */ 1124 for (i=aistart;i<aiend;i++) { 1125 row = i+r2; 1126 ierr = MatGetRow(ipmP->Ai,i,&ncols,&cols,&vals); CHKERRQ(ierr); 1127 if (ncols > 0) { 1128 /*Ai*/ 1129 ierr = MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES); CHKERRQ(ierr); 1130 /*-Ai'*/ 1131 for (j=0;j<ncols;j++) { 1132 newcol = i + c3; 1133 newrow = cols[j]; 1134 newval = -vals[j]; 1135 ierr = MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES); CHKERRQ(ierr); 1136 } 1137 } 1138 ierr = MatRestoreRow(ipmP->Ai,i,&ncols,&cols,&vals); CHKERRQ(ierr); 1139 } 1140 1141 1142 1143 /* -I */ 1144 for (i=kstart;i<kend;i++) { 1145 if (i>=r2 && i<r3) { 1146 newrow = i; 1147 newcol = i-r2+c1; 1148 newval = -1.0; 1149 MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES); CHKERRQ(ierr); 1150 } 1151 } 1152 1153 /* Copy L,Y */ 1154 ierr = VecGetOwnershipRange(ipmP->s,&sstart,&send);CHKERRQ(ierr); 1155 ierr = VecGetArray(ipmP->lamdai,&l); CHKERRQ(ierr); 1156 ierr = VecGetArray(ipmP->s,&y); CHKERRQ(ierr); 1157 1158 for (i=sstart;i<send;i++) { 1159 newcols[0] = c1+i; 1160 newcols[1] = c3+i; 1161 newvals[0] = l[i-sstart]; 1162 newvals[1] = y[i-sstart]; 1163 newrow = r3+i; 1164 ierr = MatSetValues(ipmP->K,1,&newrow,2,newcols,newvals,INSERT_VALUES); CHKERRQ(ierr); 1165 } 1166 1167 ierr = VecRestoreArray(ipmP->lamdai,&l); CHKERRQ(ierr); 1168 ierr = VecRestoreArray(ipmP->s,&y); CHKERRQ(ierr); 1169 } 1170 1171 ierr = PetscFree(indices); CHKERRQ(ierr); 1172 ierr = PetscFree(newvals); CHKERRQ(ierr); 1173 ierr = MatAssemblyBegin(ipmP->K,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1174 ierr = MatAssemblyEnd(ipmP->K,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1175 #if defined DEBUG_K 1176 PetscPrintf(PETSC_COMM_WORLD,"K\n"); MatView(ipmP->K,0); 1177 #endif 1178 PetscFunctionReturn(0); 1179 } 1180 1181 #undef __FUNCT__ 1182 #define __FUNCT__ "IPMGatherRHS" 1183 PetscErrorCode IPMGatherRHS(TaoSolver tao,Vec RHS,Vec X1,Vec X2,Vec X3,Vec X4) 1184 { 1185 TAO_IPM *ipmP = (TAO_IPM *)tao->data; 1186 PetscErrorCode ierr; 1187 PetscFunctionBegin; 1188 1189 /* rhs = [x1 (n) 1190 x2 (me) 1191 x3 (nb) 1192 x4 (nb)] */ 1193 if (X1) { 1194 ierr = VecScatterBegin(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1195 ierr = VecScatterEnd(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1196 } 1197 if (ipmP->me > 0 && X2) { 1198 ierr = VecScatterBegin(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1199 ierr = VecScatterEnd(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1200 } 1201 if (ipmP->nb > 0) { 1202 if (X3) { 1203 ierr = VecScatterBegin(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1204 ierr = VecScatterEnd(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1205 } 1206 if (X4) { 1207 ierr = VecScatterBegin(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1208 ierr = VecScatterEnd(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1209 } 1210 } 1211 #if defined(DEBUG_SCATTER) 1212 PetscPrintf(PETSC_COMM_WORLD,"X1-X4\n"); 1213 if (X1) {VecView(X1,0);} 1214 if (X2) {VecView(X2,0);} 1215 if (X3) {VecView(X3,0);} 1216 if (X4) {VecView(X4,0);} 1217 PetscPrintf(PETSC_COMM_WORLD,"RHS\n"); 1218 VecView(RHS,0); 1219 #endif 1220 PetscFunctionReturn(0); 1221 } 1222 1223 1224 1225 #undef __FUNCT__ 1226 #define __FUNCT__ "IPMScatterStep" 1227 PetscErrorCode IPMScatterStep(TaoSolver tao, Vec STEP, Vec X1, Vec X2, Vec X3, Vec X4) 1228 { 1229 TAO_IPM *ipmP = (TAO_IPM *)tao->data; 1230 PetscErrorCode ierr; 1231 PetscFunctionBegin; 1232 CHKMEMQ; 1233 /* [x1 (n) 1234 x2 (nb) may be 0 1235 x3 (me) may be 0 1236 x4 (nb) may be 0 */ 1237 if (X1) { 1238 ierr = VecScatterBegin(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1239 ierr = VecScatterEnd(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1240 } 1241 if (X2 && ipmP->nb > 0) { 1242 ierr = VecScatterBegin(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1243 ierr = VecScatterEnd(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1244 } 1245 if (X3 && ipmP->me > 0) { 1246 ierr = VecScatterBegin(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1247 ierr = VecScatterEnd(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1248 } 1249 if (X4 && ipmP->nb > 0) { 1250 ierr = VecScatterBegin(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1251 ierr = VecScatterEnd(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1252 } 1253 CHKMEMQ; 1254 #if defined(DEBUG_SCATTER) 1255 PetscPrintf(PETSC_COMM_WORLD,"Step\n"); 1256 VecView(STEP,0); 1257 PetscPrintf(PETSC_COMM_WORLD,"X1-X4\n"); 1258 if (X1) {VecView(X1,0);} 1259 if (X2) {VecView(X2,0);} 1260 if (X3) {VecView(X3,0);} 1261 if (X4) {VecView(X4,0);} 1262 #endif 1263 PetscFunctionReturn(0); 1264 } 1265 1266 1267 EXTERN_C_BEGIN 1268 1269 #undef __FUNCT__ 1270 #define __FUNCT__ "TaoCreate_IPM" 1271 PetscErrorCode TaoCreate_IPM(TaoSolver tao) 1272 { 1273 TAO_IPM *ipmP; 1274 /* const char *ipmls_type = TAOLINESEARCH_IPM; */ 1275 PetscErrorCode ierr; 1276 1277 PetscFunctionBegin; 1278 tao->ops->setup = TaoSetup_IPM; 1279 tao->ops->solve = TaoSolve_IPM; 1280 tao->ops->view = TaoView_IPM; 1281 tao->ops->setfromoptions = TaoSetFromOptions_IPM; 1282 tao->ops->destroy = TaoDestroy_IPM; 1283 /* tao->ops->computedual = TaoComputeDual_IPM; */ 1284 1285 ierr = PetscNewLog(tao,&ipmP); CHKERRQ(ierr); 1286 tao->data = (void*)ipmP; 1287 tao->max_it = 200; 1288 tao->max_funcs = 500; 1289 tao->fatol = 1e-4; 1290 tao->frtol = 1e-4; 1291 ipmP->dec = 10000; /* line search critera */ 1292 ipmP->taumin = 0.995; 1293 ipmP->monitorkkt = PETSC_FALSE; 1294 ipmP->pushs = 100; 1295 ipmP->pushnu = 100; 1296 ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp); CHKERRQ(ierr); 1297 PetscFunctionReturn(0); 1298 1299 1300 } 1301 EXTERN_C_END 1302 1303