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