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