1 #include <../src/tao/constrained/impls/ipm/pdipm.h> 2 3 /* 4 TaoPDIPMEvaluateFunctionsAndJacobians - Evaluate the objective function f, gradient fx, constraints, and all the Jacobians at current vector 5 6 Collective on tao 7 8 Input Parameter: 9 + tao - solver context 10 - x - vector at which all objects to be evaluated 11 12 Level: beginner 13 14 .seealso: `TaoPDIPMUpdateConstraints()`, `TaoPDIPMSetUpBounds()` 15 */ 16 static PetscErrorCode TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao,Vec x) 17 { 18 TAO_PDIPM *pdipm=(TAO_PDIPM*)tao->data; 19 20 PetscFunctionBegin; 21 /* Compute user objective function and gradient */ 22 PetscCall(TaoComputeObjectiveAndGradient(tao,x,&pdipm->obj,tao->gradient)); 23 24 /* Equality constraints and Jacobian */ 25 if (pdipm->Ng) { 26 PetscCall(TaoComputeEqualityConstraints(tao,x,tao->constraints_equality)); 27 PetscCall(TaoComputeJacobianEquality(tao,x,tao->jacobian_equality,tao->jacobian_equality_pre)); 28 } 29 30 /* Inequality constraints and Jacobian */ 31 if (pdipm->Nh) { 32 PetscCall(TaoComputeInequalityConstraints(tao,x,tao->constraints_inequality)); 33 PetscCall(TaoComputeJacobianInequality(tao,x,tao->jacobian_inequality,tao->jacobian_inequality_pre)); 34 } 35 PetscFunctionReturn(0); 36 } 37 38 /* 39 TaoPDIPMUpdateConstraints - Update the vectors ce and ci at x 40 41 Collective 42 43 Input Parameter: 44 + tao - Tao context 45 - x - vector at which constraints to be evaluated 46 47 Level: beginner 48 49 .seealso: `TaoPDIPMEvaluateFunctionsAndJacobians()` 50 */ 51 static PetscErrorCode TaoPDIPMUpdateConstraints(Tao tao,Vec x) 52 { 53 TAO_PDIPM *pdipm=(TAO_PDIPM*)tao->data; 54 PetscInt i,offset,offset1,k,xstart; 55 PetscScalar *carr; 56 const PetscInt *ubptr,*lbptr,*bxptr,*fxptr; 57 const PetscScalar *xarr,*xuarr,*xlarr,*garr,*harr; 58 59 PetscFunctionBegin; 60 PetscCall(VecGetOwnershipRange(x,&xstart,NULL)); 61 62 PetscCall(VecGetArrayRead(x,&xarr)); 63 PetscCall(VecGetArrayRead(tao->XU,&xuarr)); 64 PetscCall(VecGetArrayRead(tao->XL,&xlarr)); 65 66 /* (1) Update ce vector */ 67 PetscCall(VecGetArrayWrite(pdipm->ce,&carr)); 68 69 if (pdipm->Ng) { 70 /* (1.a) Inserting updated g(x) */ 71 PetscCall(VecGetArrayRead(tao->constraints_equality,&garr)); 72 PetscCall(PetscMemcpy(carr,garr,pdipm->ng*sizeof(PetscScalar))); 73 PetscCall(VecRestoreArrayRead(tao->constraints_equality,&garr)); 74 } 75 76 /* (1.b) Update xfixed */ 77 if (pdipm->Nxfixed) { 78 offset = pdipm->ng; 79 PetscCall(ISGetIndices(pdipm->isxfixed,&fxptr)); /* global indices in x */ 80 for (k=0;k < pdipm->nxfixed;k++) { 81 i = fxptr[k]-xstart; 82 carr[offset + k] = xarr[i] - xuarr[i]; 83 } 84 } 85 PetscCall(VecRestoreArrayWrite(pdipm->ce,&carr)); 86 87 /* (2) Update ci vector */ 88 PetscCall(VecGetArrayWrite(pdipm->ci,&carr)); 89 90 if (pdipm->Nh) { 91 /* (2.a) Inserting updated h(x) */ 92 PetscCall(VecGetArrayRead(tao->constraints_inequality,&harr)); 93 PetscCall(PetscMemcpy(carr,harr,pdipm->nh*sizeof(PetscScalar))); 94 PetscCall(VecRestoreArrayRead(tao->constraints_inequality,&harr)); 95 } 96 97 /* (2.b) Update xub */ 98 offset = pdipm->nh; 99 if (pdipm->Nxub) { 100 PetscCall(ISGetIndices(pdipm->isxub,&ubptr)); 101 for (k=0; k<pdipm->nxub; k++) { 102 i = ubptr[k]-xstart; 103 carr[offset + k] = xuarr[i] - xarr[i]; 104 } 105 } 106 107 if (pdipm->Nxlb) { 108 /* (2.c) Update xlb */ 109 offset += pdipm->nxub; 110 PetscCall(ISGetIndices(pdipm->isxlb,&lbptr)); /* global indices in x */ 111 for (k=0; k<pdipm->nxlb; k++) { 112 i = lbptr[k]-xstart; 113 carr[offset + k] = xarr[i] - xlarr[i]; 114 } 115 } 116 117 if (pdipm->Nxbox) { 118 /* (2.d) Update xbox */ 119 offset += pdipm->nxlb; 120 offset1 = offset + pdipm->nxbox; 121 PetscCall(ISGetIndices(pdipm->isxbox,&bxptr)); /* global indices in x */ 122 for (k=0; k<pdipm->nxbox; k++) { 123 i = bxptr[k]-xstart; /* local indices in x */ 124 carr[offset+k] = xuarr[i] - xarr[i]; 125 carr[offset1+k] = xarr[i] - xlarr[i]; 126 } 127 } 128 PetscCall(VecRestoreArrayWrite(pdipm->ci,&carr)); 129 130 /* Restoring Vectors */ 131 PetscCall(VecRestoreArrayRead(x,&xarr)); 132 PetscCall(VecRestoreArrayRead(tao->XU,&xuarr)); 133 PetscCall(VecRestoreArrayRead(tao->XL,&xlarr)); 134 PetscFunctionReturn(0); 135 } 136 137 /* 138 TaoPDIPMSetUpBounds - Create upper and lower bound vectors of x 139 140 Collective 141 142 Input Parameter: 143 . tao - holds pdipm and XL & XU 144 145 Level: beginner 146 147 .seealso: `TaoPDIPMUpdateConstraints` 148 */ 149 static PetscErrorCode TaoPDIPMSetUpBounds(Tao tao) 150 { 151 TAO_PDIPM *pdipm=(TAO_PDIPM*)tao->data; 152 const PetscScalar *xl,*xu; 153 PetscInt n,*ixlb,*ixub,*ixfixed,*ixfree,*ixbox,i,low,high,idx; 154 MPI_Comm comm; 155 PetscInt sendbuf[5],recvbuf[5]; 156 157 PetscFunctionBegin; 158 /* Creates upper and lower bounds vectors on x, if not created already */ 159 PetscCall(TaoComputeVariableBounds(tao)); 160 161 PetscCall(VecGetLocalSize(tao->XL,&n)); 162 PetscCall(PetscMalloc5(n,&ixlb,n,&ixub,n,&ixfree,n,&ixfixed,n,&ixbox)); 163 164 PetscCall(VecGetOwnershipRange(tao->XL,&low,&high)); 165 PetscCall(VecGetArrayRead(tao->XL,&xl)); 166 PetscCall(VecGetArrayRead(tao->XU,&xu)); 167 for (i=0; i<n; i++) { 168 idx = low + i; 169 if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) { 170 if (PetscRealPart(xl[i]) == PetscRealPart(xu[i])) { 171 ixfixed[pdipm->nxfixed++] = idx; 172 } else ixbox[pdipm->nxbox++] = idx; 173 } else { 174 if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) >= PETSC_INFINITY)) { 175 ixlb[pdipm->nxlb++] = idx; 176 } else if ((PetscRealPart(xl[i]) <= PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) { 177 ixub[pdipm->nxlb++] = idx; 178 } else ixfree[pdipm->nxfree++] = idx; 179 } 180 } 181 PetscCall(VecRestoreArrayRead(tao->XL,&xl)); 182 PetscCall(VecRestoreArrayRead(tao->XU,&xu)); 183 184 PetscCall(PetscObjectGetComm((PetscObject)tao,&comm)); 185 sendbuf[0] = pdipm->nxlb; 186 sendbuf[1] = pdipm->nxub; 187 sendbuf[2] = pdipm->nxfixed; 188 sendbuf[3] = pdipm->nxbox; 189 sendbuf[4] = pdipm->nxfree; 190 191 PetscCallMPI(MPI_Allreduce(sendbuf,recvbuf,5,MPIU_INT,MPI_SUM,comm)); 192 pdipm->Nxlb = recvbuf[0]; 193 pdipm->Nxub = recvbuf[1]; 194 pdipm->Nxfixed = recvbuf[2]; 195 pdipm->Nxbox = recvbuf[3]; 196 pdipm->Nxfree = recvbuf[4]; 197 198 if (pdipm->Nxlb) { 199 PetscCall(ISCreateGeneral(comm,pdipm->nxlb,ixlb,PETSC_COPY_VALUES,&pdipm->isxlb)); 200 } 201 if (pdipm->Nxub) { 202 PetscCall(ISCreateGeneral(comm,pdipm->nxub,ixub,PETSC_COPY_VALUES,&pdipm->isxub)); 203 } 204 if (pdipm->Nxfixed) { 205 PetscCall(ISCreateGeneral(comm,pdipm->nxfixed,ixfixed,PETSC_COPY_VALUES,&pdipm->isxfixed)); 206 } 207 if (pdipm->Nxbox) { 208 PetscCall(ISCreateGeneral(comm,pdipm->nxbox,ixbox,PETSC_COPY_VALUES,&pdipm->isxbox)); 209 } 210 if (pdipm->Nxfree) { 211 PetscCall(ISCreateGeneral(comm,pdipm->nxfree,ixfree,PETSC_COPY_VALUES,&pdipm->isxfree)); 212 } 213 PetscCall(PetscFree5(ixlb,ixub,ixfixed,ixbox,ixfree)); 214 PetscFunctionReturn(0); 215 } 216 217 /* 218 TaoPDIPMInitializeSolution - Initialize PDIPM solution X = [x; lambdae; lambdai; z]. 219 X consists of four subvectors in the order [x; lambdae; lambdai; z]. These 220 four subvectors need to be initialized and its values copied over to X. Instead 221 of copying, we use VecPlace/ResetArray functions to share the memory locations for 222 X and the subvectors 223 224 Collective 225 226 Input Parameter: 227 . tao - Tao context 228 229 Level: beginner 230 */ 231 static PetscErrorCode TaoPDIPMInitializeSolution(Tao tao) 232 { 233 TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 234 PetscScalar *Xarr,*z,*lambdai; 235 PetscInt i; 236 const PetscScalar *xarr,*h; 237 238 PetscFunctionBegin; 239 PetscCall(VecGetArrayWrite(pdipm->X,&Xarr)); 240 241 /* Set Initialize X.x = tao->solution */ 242 PetscCall(VecGetArrayRead(tao->solution,&xarr)); 243 PetscCall(PetscMemcpy(Xarr,xarr,pdipm->nx*sizeof(PetscScalar))); 244 PetscCall(VecRestoreArrayRead(tao->solution,&xarr)); 245 246 /* Initialize X.lambdae = 0.0 */ 247 if (pdipm->lambdae) PetscCall(VecSet(pdipm->lambdae,0.0)); 248 249 /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */ 250 if (pdipm->Nci) { 251 PetscCall(VecSet(pdipm->lambdai,pdipm->push_init_lambdai)); 252 PetscCall(VecSet(pdipm->z,pdipm->push_init_slack)); 253 254 /* Additional modification for X.lambdai and X.z */ 255 PetscCall(VecGetArrayWrite(pdipm->lambdai,&lambdai)); 256 PetscCall(VecGetArrayWrite(pdipm->z,&z)); 257 if (pdipm->Nh) { 258 PetscCall(VecGetArrayRead(tao->constraints_inequality,&h)); 259 for (i=0; i < pdipm->nh; i++) { 260 if (h[i] < -pdipm->push_init_slack) z[i] = -h[i]; 261 if (pdipm->mu/z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu/z[i]; 262 } 263 PetscCall(VecRestoreArrayRead(tao->constraints_inequality,&h)); 264 } 265 PetscCall(VecRestoreArrayWrite(pdipm->lambdai,&lambdai)); 266 PetscCall(VecRestoreArrayWrite(pdipm->z,&z)); 267 } 268 269 PetscCall(VecRestoreArrayWrite(pdipm->X,&Xarr)); 270 PetscFunctionReturn(0); 271 } 272 273 /* 274 TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X 275 276 Input Parameter: 277 snes - SNES context 278 X - KKT Vector 279 *ctx - pdipm context 280 281 Output Parameter: 282 J - Hessian matrix 283 Jpre - Preconditioner 284 */ 285 static PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes,Vec X, Mat J, Mat Jpre, void *ctx) 286 { 287 Tao tao=(Tao)ctx; 288 TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 289 PetscInt i,row,cols[2],Jrstart,rjstart,nc,j; 290 const PetscInt *aj,*ranges,*Jranges,*rranges,*cranges; 291 const PetscScalar *Xarr,*aa; 292 PetscScalar vals[2]; 293 PetscInt proc,nx_all,*nce_all=pdipm->nce_all; 294 MPI_Comm comm; 295 PetscMPIInt rank,size; 296 297 PetscFunctionBegin; 298 PetscCall(PetscObjectGetComm((PetscObject)snes,&comm)); 299 PetscCallMPI(MPI_Comm_rank(comm,&rank)); 300 PetscCallMPI(MPI_Comm_rank(comm,&size)); 301 302 PetscCall(MatGetOwnershipRanges(Jpre,&Jranges)); 303 PetscCall(MatGetOwnershipRange(Jpre,&Jrstart,NULL)); 304 PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&rranges)); 305 PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&cranges)); 306 307 PetscCall(VecGetArrayRead(X,&Xarr)); 308 309 /* (1) insert Z and Ci to the 4th block of Jpre -- overwrite existing values */ 310 if (pdipm->solve_symmetric_kkt) { /* 1 for eq 17 revised pdipm doc 0 for eq 18 (symmetric KKT) */ 311 vals[0] = 1.0; 312 for (i=0; i < pdipm->nci; i++) { 313 row = Jrstart + pdipm->off_z + i; 314 cols[0] = Jrstart + pdipm->off_lambdai + i; 315 cols[1] = row; 316 vals[1] = Xarr[pdipm->off_lambdai + i]/Xarr[pdipm->off_z + i]; 317 PetscCall(MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES)); 318 } 319 } else { 320 for (i=0; i < pdipm->nci; i++) { 321 row = Jrstart + pdipm->off_z + i; 322 cols[0] = Jrstart + pdipm->off_lambdai + i; 323 cols[1] = row; 324 vals[0] = Xarr[pdipm->off_z + i]; 325 vals[1] = Xarr[pdipm->off_lambdai + i]; 326 PetscCall(MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES)); 327 } 328 } 329 330 /* (2) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */ 331 if (pdipm->Ng) { 332 PetscCall(MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL)); 333 for (i=0; i<pdipm->ng; i++) { 334 row = Jrstart + pdipm->off_lambdae + i; 335 336 PetscCall(MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa)); 337 proc = 0; 338 for (j=0; j < nc; j++) { 339 while (aj[j] >= cranges[proc+1]) proc++; 340 cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 341 PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 342 } 343 PetscCall(MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa)); 344 if (pdipm->kkt_pd) { 345 /* add shift \delta_c */ 346 PetscCall(MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES)); 347 } 348 } 349 } 350 351 /* (3) insert 3rd row block of Jpre: [ -grad h, 0, deltac, I] */ 352 if (pdipm->Nh) { 353 PetscCall(MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL)); 354 for (i=0; i < pdipm->nh; i++) { 355 row = Jrstart + pdipm->off_lambdai + i; 356 PetscCall(MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa)); 357 proc = 0; 358 for (j=0; j < nc; j++) { 359 while (aj[j] >= cranges[proc+1]) proc++; 360 cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 361 PetscCall(MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES)); 362 } 363 PetscCall(MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa)); 364 if (pdipm->kkt_pd) { 365 /* add shift \delta_c */ 366 PetscCall(MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES)); 367 } 368 } 369 } 370 371 /* (4) insert 1st row block of Jpre: [Wxx, grad g', -grad h', 0] */ 372 if (pdipm->Ng) { /* grad g' */ 373 PetscCall(MatTranspose(tao->jacobian_equality,MAT_REUSE_MATRIX,&pdipm->jac_equality_trans)); 374 } 375 if (pdipm->Nh) { /* grad h' */ 376 PetscCall(MatTranspose(tao->jacobian_inequality,MAT_REUSE_MATRIX,&pdipm->jac_inequality_trans)); 377 } 378 379 PetscCall(VecPlaceArray(pdipm->x,Xarr)); 380 PetscCall(TaoComputeHessian(tao,pdipm->x,tao->hessian,tao->hessian_pre)); 381 PetscCall(VecResetArray(pdipm->x)); 382 383 PetscCall(MatGetOwnershipRange(tao->hessian,&rjstart,NULL)); 384 for (i=0; i<pdipm->nx; i++) { 385 row = Jrstart + i; 386 387 /* insert Wxx = fxx + ... -- provided by user */ 388 PetscCall(MatGetRow(tao->hessian,i+rjstart,&nc,&aj,&aa)); 389 proc = 0; 390 for (j=0; j < nc; j++) { 391 while (aj[j] >= cranges[proc+1]) proc++; 392 cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 393 if (row == cols[0] && pdipm->kkt_pd) { 394 /* add shift deltaw to Wxx component */ 395 PetscCall(MatSetValue(Jpre,row,cols[0],aa[j]+pdipm->deltaw,INSERT_VALUES)); 396 } else { 397 PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 398 } 399 } 400 PetscCall(MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,&aa)); 401 402 /* insert grad g' */ 403 if (pdipm->ng) { 404 PetscCall(MatGetRow(pdipm->jac_equality_trans,i+rjstart,&nc,&aj,&aa)); 405 PetscCall(MatGetOwnershipRanges(tao->jacobian_equality,&ranges)); 406 proc = 0; 407 for (j=0; j < nc; j++) { 408 /* find row ownership of */ 409 while (aj[j] >= ranges[proc+1]) proc++; 410 nx_all = rranges[proc+1] - rranges[proc]; 411 cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 412 PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 413 } 414 PetscCall(MatRestoreRow(pdipm->jac_equality_trans,i+rjstart,&nc,&aj,&aa)); 415 } 416 417 /* insert -grad h' */ 418 if (pdipm->nh) { 419 PetscCall(MatGetRow(pdipm->jac_inequality_trans,i+rjstart,&nc,&aj,&aa)); 420 PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality,&ranges)); 421 proc = 0; 422 for (j=0; j < nc; j++) { 423 /* find row ownership of */ 424 while (aj[j] >= ranges[proc+1]) proc++; 425 nx_all = rranges[proc+1] - rranges[proc]; 426 cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 427 PetscCall(MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES)); 428 } 429 PetscCall(MatRestoreRow(pdipm->jac_inequality_trans,i+rjstart,&nc,&aj,&aa)); 430 } 431 } 432 PetscCall(VecRestoreArrayRead(X,&Xarr)); 433 434 /* (6) assemble Jpre and J */ 435 PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 436 PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 437 438 if (J != Jpre) { 439 PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 440 PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 441 } 442 PetscFunctionReturn(0); 443 } 444 445 /* 446 TaoSnesFunction_PDIPM - Evaluate KKT function at X 447 448 Input Parameter: 449 snes - SNES context 450 X - KKT Vector 451 *ctx - pdipm 452 453 Output Parameter: 454 F - Updated Lagrangian vector 455 */ 456 static PetscErrorCode TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,void *ctx) 457 { 458 Tao tao=(Tao)ctx; 459 TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 460 PetscScalar *Farr; 461 Vec x,L1; 462 PetscInt i; 463 const PetscScalar *Xarr,*carr,*zarr,*larr; 464 465 PetscFunctionBegin; 466 PetscCall(VecSet(F,0.0)); 467 468 PetscCall(VecGetArrayRead(X,&Xarr)); 469 PetscCall(VecGetArrayWrite(F,&Farr)); 470 471 /* (0) Evaluate f, fx, gradG, gradH at X.x Note: pdipm->x is not changed below */ 472 x = pdipm->x; 473 PetscCall(VecPlaceArray(x,Xarr)); 474 PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao,x)); 475 476 /* Update ce, ci, and Jci at X.x */ 477 PetscCall(TaoPDIPMUpdateConstraints(tao,x)); 478 PetscCall(VecResetArray(x)); 479 480 /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */ 481 L1 = pdipm->x; 482 PetscCall(VecPlaceArray(L1,Farr)); /* L1 = 0.0 */ 483 if (pdipm->Nci) { 484 if (pdipm->Nh) { 485 /* L1 += gradH'*DI. Note: tao->DI is not changed below */ 486 PetscCall(VecPlaceArray(tao->DI,Xarr+pdipm->off_lambdai)); 487 PetscCall(MatMultTransposeAdd(tao->jacobian_inequality,tao->DI,L1,L1)); 488 PetscCall(VecResetArray(tao->DI)); 489 } 490 491 /* L1 += Jci_xb'*lambdai_xb */ 492 PetscCall(VecPlaceArray(pdipm->lambdai_xb,Xarr+pdipm->off_lambdai+pdipm->nh)); 493 PetscCall(MatMultTransposeAdd(pdipm->Jci_xb,pdipm->lambdai_xb,L1,L1)); 494 PetscCall(VecResetArray(pdipm->lambdai_xb)); 495 496 /* L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */ 497 PetscCall(VecScale(L1,-1.0)); 498 } 499 500 /* L1 += fx */ 501 PetscCall(VecAXPY(L1,1.0,tao->gradient)); 502 503 if (pdipm->Nce) { 504 if (pdipm->Ng) { 505 /* L1 += gradG'*DE. Note: tao->DE is not changed below */ 506 PetscCall(VecPlaceArray(tao->DE,Xarr+pdipm->off_lambdae)); 507 PetscCall(MatMultTransposeAdd(tao->jacobian_equality,tao->DE,L1,L1)); 508 PetscCall(VecResetArray(tao->DE)); 509 } 510 if (pdipm->Nxfixed) { 511 /* L1 += Jce_xfixed'*lambdae_xfixed */ 512 PetscCall(VecPlaceArray(pdipm->lambdae_xfixed,Xarr+pdipm->off_lambdae+pdipm->ng)); 513 PetscCall(MatMultTransposeAdd(pdipm->Jce_xfixed,pdipm->lambdae_xfixed,L1,L1)); 514 PetscCall(VecResetArray(pdipm->lambdae_xfixed)); 515 } 516 } 517 PetscCall(VecResetArray(L1)); 518 519 /* (2) L2 = ce(x) */ 520 if (pdipm->Nce) { 521 PetscCall(VecGetArrayRead(pdipm->ce,&carr)); 522 for (i=0; i<pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i]; 523 PetscCall(VecRestoreArrayRead(pdipm->ce,&carr)); 524 } 525 526 if (pdipm->Nci) { 527 if (pdipm->solve_symmetric_kkt) { 528 /* (3) L3 = z - ci(x); 529 (4) L4 = Lambdai * e - mu/z *e */ 530 PetscCall(VecGetArrayRead(pdipm->ci,&carr)); 531 larr = Xarr+pdipm->off_lambdai; 532 zarr = Xarr+pdipm->off_z; 533 for (i=0; i<pdipm->nci; i++) { 534 Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 535 Farr[pdipm->off_z + i] = larr[i] - pdipm->mu/zarr[i]; 536 } 537 PetscCall(VecRestoreArrayRead(pdipm->ci,&carr)); 538 } else { 539 /* (3) L3 = z - ci(x); 540 (4) L4 = Z * Lambdai * e - mu * e */ 541 PetscCall(VecGetArrayRead(pdipm->ci,&carr)); 542 larr = Xarr+pdipm->off_lambdai; 543 zarr = Xarr+pdipm->off_z; 544 for (i=0; i<pdipm->nci; i++) { 545 Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 546 Farr[pdipm->off_z + i] = zarr[i]*larr[i] - pdipm->mu; 547 } 548 PetscCall(VecRestoreArrayRead(pdipm->ci,&carr)); 549 } 550 } 551 552 PetscCall(VecRestoreArrayRead(X,&Xarr)); 553 PetscCall(VecRestoreArrayWrite(F,&Farr)); 554 PetscFunctionReturn(0); 555 } 556 557 /* 558 Evaluate F(X); then update update tao->gnorm0, tao->step = mu, 559 tao->residual = norm2(F_x,F_z) and tao->cnorm = norm2(F_ce,F_ci). 560 */ 561 static PetscErrorCode TaoSNESFunction_PDIPM_residual(SNES snes,Vec X,Vec F,void *ctx) 562 { 563 Tao tao=(Tao)ctx; 564 TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 565 PetscScalar *Farr,*tmparr; 566 Vec L1; 567 PetscInt i; 568 PetscReal res[2],cnorm[2]; 569 const PetscScalar *Xarr=NULL; 570 571 PetscFunctionBegin; 572 PetscCall(TaoSNESFunction_PDIPM(snes,X,F,(void*)tao)); 573 PetscCall(VecGetArrayWrite(F,&Farr)); 574 PetscCall(VecGetArrayRead(X,&Xarr)); 575 576 /* compute res[0] = norm2(F_x) */ 577 L1 = pdipm->x; 578 PetscCall(VecPlaceArray(L1,Farr)); 579 PetscCall(VecNorm(L1,NORM_2,&res[0])); 580 PetscCall(VecResetArray(L1)); 581 582 /* compute res[1] = norm2(F_z), cnorm[1] = norm2(F_ci) */ 583 if (pdipm->z) { 584 if (pdipm->solve_symmetric_kkt) { 585 PetscCall(VecPlaceArray(pdipm->z,Farr+pdipm->off_z)); 586 if (pdipm->Nci) { 587 PetscCall(VecGetArrayWrite(pdipm->z,&tmparr)); 588 for (i=0; i<pdipm->nci; i++) tmparr[i] *= Xarr[pdipm->off_z + i]; 589 PetscCall(VecRestoreArrayWrite(pdipm->z,&tmparr)); 590 } 591 592 PetscCall(VecNorm(pdipm->z,NORM_2,&res[1])); 593 594 if (pdipm->Nci) { 595 PetscCall(VecGetArrayWrite(pdipm->z,&tmparr)); 596 for (i=0; i<pdipm->nci; i++) { 597 tmparr[i] /= Xarr[pdipm->off_z + i]; 598 } 599 PetscCall(VecRestoreArrayWrite(pdipm->z,&tmparr)); 600 } 601 PetscCall(VecResetArray(pdipm->z)); 602 } else { /* !solve_symmetric_kkt */ 603 PetscCall(VecPlaceArray(pdipm->z,Farr+pdipm->off_z)); 604 PetscCall(VecNorm(pdipm->z,NORM_2,&res[1])); 605 PetscCall(VecResetArray(pdipm->z)); 606 } 607 608 PetscCall(VecPlaceArray(pdipm->ci,Farr+pdipm->off_lambdai)); 609 PetscCall(VecNorm(pdipm->ci,NORM_2,&cnorm[1])); 610 PetscCall(VecResetArray(pdipm->ci)); 611 } else { 612 res[1] = 0.0; cnorm[1] = 0.0; 613 } 614 615 /* compute cnorm[0] = norm2(F_ce) */ 616 if (pdipm->Nce) { 617 PetscCall(VecPlaceArray(pdipm->ce,Farr+pdipm->off_lambdae)); 618 PetscCall(VecNorm(pdipm->ce,NORM_2,&cnorm[0])); 619 PetscCall(VecResetArray(pdipm->ce)); 620 } else cnorm[0] = 0.0; 621 622 PetscCall(VecRestoreArrayWrite(F,&Farr)); 623 PetscCall(VecRestoreArrayRead(X,&Xarr)); 624 625 tao->gnorm0 = tao->residual; 626 tao->residual = PetscSqrtReal(res[0]*res[0] + res[1]*res[1]); 627 tao->cnorm = PetscSqrtReal(cnorm[0]*cnorm[0] + cnorm[1]*cnorm[1]); 628 tao->step = pdipm->mu; 629 PetscFunctionReturn(0); 630 } 631 632 /* 633 KKTAddShifts - Check the inertia of Cholesky factor of KKT matrix. 634 If it does not match the numbers of prime and dual variables, add shifts to the KKT matrix. 635 */ 636 static PetscErrorCode KKTAddShifts(Tao tao,SNES snes,Vec X) 637 { 638 TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 639 KSP ksp; 640 PC pc; 641 Mat Factor; 642 PetscBool isCHOL; 643 PetscInt nneg,nzero,npos; 644 645 PetscFunctionBegin; 646 /* Get the inertia of Cholesky factor */ 647 PetscCall(SNESGetKSP(snes,&ksp)); 648 PetscCall(KSPGetPC(ksp,&pc)); 649 PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL)); 650 if (!isCHOL) PetscFunctionReturn(0); 651 652 PetscCall(PCFactorGetMatrix(pc,&Factor)); 653 PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 654 655 if (npos < pdipm->Nx+pdipm->Nci) { 656 pdipm->deltaw = PetscMax(pdipm->lastdeltaw/3, 1.e-4*PETSC_MACHINE_EPSILON); 657 PetscCall(PetscInfo(tao,"Test reduced deltaw=%g; previous MatInertia: nneg %" PetscInt_FMT ", nzero %" PetscInt_FMT ", npos %" PetscInt_FMT "(<%" PetscInt_FMT ")\n",(double)pdipm->deltaw,nneg,nzero,npos,pdipm->Nx+pdipm->Nci)); 658 PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 659 PetscCall(PCSetUp(pc)); 660 PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 661 662 if (npos < pdipm->Nx+pdipm->Nci) { 663 pdipm->deltaw = pdipm->lastdeltaw; /* in case reduction update does not help, this prevents that step from impacting increasing update */ 664 while (npos < pdipm->Nx+pdipm->Nci && pdipm->deltaw <= 1./PETSC_SMALL) { /* increase deltaw */ 665 PetscCall(PetscInfo(tao," deltaw=%g fails, MatInertia: nneg %" PetscInt_FMT ", nzero %" PetscInt_FMT ", npos %" PetscInt_FMT "(<%" PetscInt_FMT ")\n",(double)pdipm->deltaw,nneg,nzero,npos,pdipm->Nx+pdipm->Nci)); 666 pdipm->deltaw = PetscMin(8*pdipm->deltaw,PetscPowReal(10,20)); 667 PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 668 PetscCall(PCSetUp(pc)); 669 PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 670 } 671 672 PetscCheck(pdipm->deltaw < 1./PETSC_SMALL,PetscObjectComm((PetscObject)tao),PETSC_ERR_CONV_FAILED,"Reached maximum delta w will not converge, try different initial x0"); 673 674 PetscCall(PetscInfo(tao,"Updated deltaw %g\n",(double)pdipm->deltaw)); 675 pdipm->lastdeltaw = pdipm->deltaw; 676 pdipm->deltaw = 0.0; 677 } 678 } 679 680 if (nzero) { /* Jacobian is singular */ 681 if (pdipm->deltac == 0.0) { 682 pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON; 683 } else { 684 pdipm->deltac = pdipm->deltac*PetscPowReal(pdipm->mu,.25); 685 } 686 PetscCall(PetscInfo(tao,"Updated deltac=%g, MatInertia: nneg %" PetscInt_FMT ", nzero %" PetscInt_FMT "(!=0), npos %" PetscInt_FMT "\n",(double)pdipm->deltac,nneg,nzero,npos)); 687 PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 688 PetscCall(PCSetUp(pc)); 689 PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 690 } 691 PetscFunctionReturn(0); 692 } 693 694 /* 695 PCPreSolve_PDIPM -- called betwee MatFactorNumeric() and MatSolve() 696 */ 697 PetscErrorCode PCPreSolve_PDIPM(PC pc,KSP ksp) 698 { 699 Tao tao; 700 TAO_PDIPM *pdipm; 701 702 PetscFunctionBegin; 703 PetscCall(KSPGetApplicationContext(ksp,&tao)); 704 pdipm = (TAO_PDIPM*)tao->data; 705 PetscCall(KKTAddShifts(tao,pdipm->snes,pdipm->X)); 706 PetscFunctionReturn(0); 707 } 708 709 /* 710 SNESLineSearch_PDIPM - Custom line search used with PDIPM. 711 712 Collective on TAO 713 714 Notes: 715 This routine employs a simple backtracking line-search to keep 716 the slack variables (z) and inequality constraints Lagrange multipliers 717 (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars 718 alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the 719 slack variables are updated as X = X - alpha_d*dx. The constraint multipliers 720 are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu 721 is also updated as mu = mu + z'lambdai/Nci 722 */ 723 static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch,void *ctx) 724 { 725 Tao tao=(Tao)ctx; 726 TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 727 SNES snes; 728 Vec X,F,Y; 729 PetscInt i,iter; 730 PetscReal alpha_p=1.0,alpha_d=1.0,alpha[4]; 731 PetscScalar *Xarr,*z,*lambdai,dot,*taosolarr; 732 const PetscScalar *dXarr,*dz,*dlambdai; 733 734 PetscFunctionBegin; 735 PetscCall(SNESLineSearchGetSNES(linesearch,&snes)); 736 PetscCall(SNESGetIterationNumber(snes,&iter)); 737 738 PetscCall(SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_SUCCEEDED)); 739 PetscCall(SNESLineSearchGetVecs(linesearch,&X,&F,&Y,NULL,NULL)); 740 741 PetscCall(VecGetArrayWrite(X,&Xarr)); 742 PetscCall(VecGetArrayRead(Y,&dXarr)); 743 z = Xarr + pdipm->off_z; 744 dz = dXarr + pdipm->off_z; 745 for (i=0; i < pdipm->nci; i++) { 746 if (z[i] - dz[i] < 0.0) alpha_p = PetscMin(alpha_p, 0.9999*z[i]/dz[i]); 747 } 748 749 lambdai = Xarr + pdipm->off_lambdai; 750 dlambdai = dXarr + pdipm->off_lambdai; 751 752 for (i=0; i<pdipm->nci; i++) { 753 if (lambdai[i] - dlambdai[i] < 0.0) alpha_d = PetscMin(0.9999*lambdai[i]/dlambdai[i], alpha_d); 754 } 755 756 alpha[0] = alpha_p; 757 alpha[1] = alpha_d; 758 PetscCall(VecRestoreArrayRead(Y,&dXarr)); 759 PetscCall(VecRestoreArrayWrite(X,&Xarr)); 760 761 /* alpha = min(alpha) over all processes */ 762 PetscCallMPI(MPI_Allreduce(alpha,alpha+2,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)tao))); 763 764 alpha_p = alpha[2]; 765 alpha_d = alpha[3]; 766 767 /* X = X - alpha * Y */ 768 PetscCall(VecGetArrayWrite(X,&Xarr)); 769 PetscCall(VecGetArrayRead(Y,&dXarr)); 770 for (i=0; i<pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i]; 771 for (i=0; i<pdipm->nce; i++) Xarr[i+pdipm->off_lambdae] -= alpha_d * dXarr[i+pdipm->off_lambdae]; 772 773 for (i=0; i<pdipm->nci; i++) { 774 Xarr[i+pdipm->off_lambdai] -= alpha_d * dXarr[i+pdipm->off_lambdai]; 775 Xarr[i+pdipm->off_z] -= alpha_p * dXarr[i+pdipm->off_z]; 776 } 777 PetscCall(VecGetArrayWrite(tao->solution,&taosolarr)); 778 PetscCall(PetscMemcpy(taosolarr,Xarr,pdipm->nx*sizeof(PetscScalar))); 779 PetscCall(VecRestoreArrayWrite(tao->solution,&taosolarr)); 780 781 PetscCall(VecRestoreArrayWrite(X,&Xarr)); 782 PetscCall(VecRestoreArrayRead(Y,&dXarr)); 783 784 /* Update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */ 785 if (pdipm->z) { 786 PetscCall(VecDot(pdipm->z,pdipm->lambdai,&dot)); 787 } else dot = 0.0; 788 789 /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu) */ 790 pdipm->mu = pdipm->mu_update_factor * dot/pdipm->Nci; 791 792 /* Update F; get tao->residual and tao->cnorm */ 793 PetscCall(TaoSNESFunction_PDIPM_residual(snes,X,F,(void*)tao)); 794 795 tao->niter++; 796 PetscCall(TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter)); 797 PetscCall(TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu)); 798 799 PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 800 if (tao->reason) PetscCall(SNESSetConvergedReason(snes,SNES_CONVERGED_FNORM_ABS)); 801 PetscFunctionReturn(0); 802 } 803 804 /* 805 TaoSolve_PDIPM 806 807 Input Parameter: 808 tao - TAO context 809 810 Output Parameter: 811 tao - TAO context 812 */ 813 PetscErrorCode TaoSolve_PDIPM(Tao tao) 814 { 815 TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 816 SNESLineSearch linesearch; /* SNESLineSearch context */ 817 Vec dummy; 818 819 PetscFunctionBegin; 820 PetscCheck(tao->constraints_equality || tao->constraints_inequality,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_NULL,"Equality and inequality constraints are not set. Either set them or switch to a different algorithm"); 821 822 /* Initialize all variables */ 823 PetscCall(TaoPDIPMInitializeSolution(tao)); 824 825 /* Set linesearch */ 826 PetscCall(SNESGetLineSearch(pdipm->snes,&linesearch)); 827 PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHSHELL)); 828 PetscCall(SNESLineSearchShellSetUserFunc(linesearch,SNESLineSearch_PDIPM,tao)); 829 PetscCall(SNESLineSearchSetFromOptions(linesearch)); 830 831 tao->reason = TAO_CONTINUE_ITERATING; 832 833 /* -tao_monitor for iteration 0 and check convergence */ 834 PetscCall(VecDuplicate(pdipm->X,&dummy)); 835 PetscCall(TaoSNESFunction_PDIPM_residual(pdipm->snes,pdipm->X,dummy,(void*)tao)); 836 837 PetscCall(TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter)); 838 PetscCall(TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu)); 839 PetscCall(VecDestroy(&dummy)); 840 PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 841 if (tao->reason) PetscCall(SNESSetConvergedReason(pdipm->snes,SNES_CONVERGED_FNORM_ABS)); 842 843 while (tao->reason == TAO_CONTINUE_ITERATING) { 844 SNESConvergedReason reason; 845 PetscCall(SNESSolve(pdipm->snes,NULL,pdipm->X)); 846 847 /* Check SNES convergence */ 848 PetscCall(SNESGetConvergedReason(pdipm->snes,&reason)); 849 if (reason < 0) { 850 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes),"SNES solve did not converged due to reason %s\n",SNESConvergedReasons[reason])); 851 } 852 853 /* Check TAO convergence */ 854 PetscCheck(!PetscIsInfOrNanReal(pdipm->obj),PETSC_COMM_SELF,PETSC_ERR_SUP,"User-provided compute function generated Inf or NaN"); 855 } 856 PetscFunctionReturn(0); 857 } 858 859 /* 860 TaoView_PDIPM - View PDIPM 861 862 Input Parameter: 863 tao - TAO object 864 viewer - PetscViewer 865 866 Output: 867 */ 868 PetscErrorCode TaoView_PDIPM(Tao tao,PetscViewer viewer) 869 { 870 TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 871 872 PetscFunctionBegin; 873 tao->constrained = PETSC_TRUE; 874 PetscCall(PetscViewerASCIIPushTab(viewer)); 875 PetscCall(PetscViewerASCIIPrintf(viewer,"Number of prime=%" PetscInt_FMT ", Number of dual=%" PetscInt_FMT "\n",pdipm->Nx+pdipm->Nci,pdipm->Nce + pdipm->Nci)); 876 if (pdipm->kkt_pd) { 877 PetscCall(PetscViewerASCIIPrintf(viewer,"KKT shifts deltaw=%g, deltac=%g\n",(double)pdipm->deltaw,(double)pdipm->deltac)); 878 } 879 PetscCall(PetscViewerASCIIPopTab(viewer)); 880 PetscFunctionReturn(0); 881 } 882 883 /* 884 TaoSetup_PDIPM - Sets up tao and pdipm 885 886 Input Parameter: 887 tao - TAO object 888 889 Output: pdipm - initialized object 890 */ 891 PetscErrorCode TaoSetup_PDIPM(Tao tao) 892 { 893 TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 894 MPI_Comm comm; 895 PetscMPIInt size; 896 PetscInt row,col,Jcrstart,Jcrend,k,tmp,nc,proc,*nh_all,*ng_all; 897 PetscInt offset,*xa,*xb,i,j,rstart,rend; 898 PetscScalar one=1.0,neg_one=-1.0; 899 const PetscInt *cols,*rranges,*cranges,*aj,*ranges; 900 const PetscScalar *aa,*Xarr; 901 Mat J; 902 Mat Jce_xfixed_trans,Jci_xb_trans; 903 PetscInt *dnz,*onz,rjstart,nx_all,*nce_all,*Jranges,cols1[2]; 904 905 PetscFunctionBegin; 906 PetscCall(PetscObjectGetComm((PetscObject)tao,&comm)); 907 PetscCallMPI(MPI_Comm_size(comm,&size)); 908 909 /* (1) Setup Bounds and create Tao vectors */ 910 PetscCall(TaoPDIPMSetUpBounds(tao)); 911 912 if (!tao->gradient) { 913 PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 914 PetscCall(VecDuplicate(tao->solution,&tao->stepdirection)); 915 } 916 917 /* (2) Get sizes */ 918 /* Size of vector x - This is set by TaoSetSolution */ 919 PetscCall(VecGetSize(tao->solution,&pdipm->Nx)); 920 PetscCall(VecGetLocalSize(tao->solution,&pdipm->nx)); 921 922 /* Size of equality constraints and vectors */ 923 if (tao->constraints_equality) { 924 PetscCall(VecGetSize(tao->constraints_equality,&pdipm->Ng)); 925 PetscCall(VecGetLocalSize(tao->constraints_equality,&pdipm->ng)); 926 } else { 927 pdipm->ng = pdipm->Ng = 0; 928 } 929 930 pdipm->nce = pdipm->ng + pdipm->nxfixed; 931 pdipm->Nce = pdipm->Ng + pdipm->Nxfixed; 932 933 /* Size of inequality constraints and vectors */ 934 if (tao->constraints_inequality) { 935 PetscCall(VecGetSize(tao->constraints_inequality,&pdipm->Nh)); 936 PetscCall(VecGetLocalSize(tao->constraints_inequality,&pdipm->nh)); 937 } else { 938 pdipm->nh = pdipm->Nh = 0; 939 } 940 941 pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2*pdipm->nxbox; 942 pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2*pdipm->Nxbox; 943 944 /* Full size of the KKT system to be solved */ 945 pdipm->n = pdipm->nx + pdipm->nce + 2*pdipm->nci; 946 pdipm->N = pdipm->Nx + pdipm->Nce + 2*pdipm->Nci; 947 948 /* (3) Offsets for subvectors */ 949 pdipm->off_lambdae = pdipm->nx; 950 pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce; 951 pdipm->off_z = pdipm->off_lambdai + pdipm->nci; 952 953 /* (4) Create vectors and subvectors */ 954 /* Ce and Ci vectors */ 955 PetscCall(VecCreate(comm,&pdipm->ce)); 956 PetscCall(VecSetSizes(pdipm->ce,pdipm->nce,pdipm->Nce)); 957 PetscCall(VecSetFromOptions(pdipm->ce)); 958 959 PetscCall(VecCreate(comm,&pdipm->ci)); 960 PetscCall(VecSetSizes(pdipm->ci,pdipm->nci,pdipm->Nci)); 961 PetscCall(VecSetFromOptions(pdipm->ci)); 962 963 /* X=[x; lambdae; lambdai; z] for the big KKT system */ 964 PetscCall(VecCreate(comm,&pdipm->X)); 965 PetscCall(VecSetSizes(pdipm->X,pdipm->n,pdipm->N)); 966 PetscCall(VecSetFromOptions(pdipm->X)); 967 968 /* Subvectors; they share local arrays with X */ 969 PetscCall(VecGetArrayRead(pdipm->X,&Xarr)); 970 /* x shares local array with X.x */ 971 if (pdipm->Nx) { 972 PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nx,pdipm->Nx,Xarr,&pdipm->x)); 973 } 974 975 /* lambdae shares local array with X.lambdae */ 976 if (pdipm->Nce) { 977 PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nce,pdipm->Nce,Xarr+pdipm->off_lambdae,&pdipm->lambdae)); 978 } 979 980 /* tao->DE shares local array with X.lambdae_g */ 981 if (pdipm->Ng) { 982 PetscCall(VecCreateMPIWithArray(comm,1,pdipm->ng,pdipm->Ng,Xarr+pdipm->off_lambdae,&tao->DE)); 983 984 PetscCall(VecCreate(comm,&pdipm->lambdae_xfixed)); 985 PetscCall(VecSetSizes(pdipm->lambdae_xfixed,pdipm->nxfixed,PETSC_DECIDE)); 986 PetscCall(VecSetFromOptions(pdipm->lambdae_xfixed)); 987 } 988 989 if (pdipm->Nci) { 990 /* lambdai shares local array with X.lambdai */ 991 PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_lambdai,&pdipm->lambdai)); 992 993 /* z for slack variables; it shares local array with X.z */ 994 PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_z,&pdipm->z)); 995 } 996 997 /* tao->DI which shares local array with X.lambdai_h */ 998 if (pdipm->Nh) { 999 PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nh,pdipm->Nh,Xarr+pdipm->off_lambdai,&tao->DI)); 1000 } 1001 PetscCall(VecCreate(comm,&pdipm->lambdai_xb)); 1002 PetscCall(VecSetSizes(pdipm->lambdai_xb,(pdipm->nci - pdipm->nh),PETSC_DECIDE)); 1003 PetscCall(VecSetFromOptions(pdipm->lambdai_xb)); 1004 1005 PetscCall(VecRestoreArrayRead(pdipm->X,&Xarr)); 1006 1007 /* (5) Create Jacobians Jce_xfixed and Jci */ 1008 /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */ 1009 if (pdipm->Nxfixed) { 1010 /* Create Jce_xfixed */ 1011 PetscCall(MatCreate(comm,&pdipm->Jce_xfixed)); 1012 PetscCall(MatSetSizes(pdipm->Jce_xfixed,pdipm->nxfixed,pdipm->nx,PETSC_DECIDE,pdipm->Nx)); 1013 PetscCall(MatSetFromOptions(pdipm->Jce_xfixed)); 1014 PetscCall(MatSeqAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL)); 1015 PetscCall(MatMPIAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL,1,NULL)); 1016 1017 PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,&Jcrend)); 1018 PetscCall(ISGetIndices(pdipm->isxfixed,&cols)); 1019 k = 0; 1020 for (row = Jcrstart; row < Jcrend; row++) { 1021 PetscCall(MatSetValues(pdipm->Jce_xfixed,1,&row,1,cols+k,&one,INSERT_VALUES)); 1022 k++; 1023 } 1024 PetscCall(ISRestoreIndices(pdipm->isxfixed, &cols)); 1025 PetscCall(MatAssemblyBegin(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY)); 1026 PetscCall(MatAssemblyEnd(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY)); 1027 } 1028 1029 /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */ 1030 PetscCall(MatCreate(comm,&pdipm->Jci_xb)); 1031 PetscCall(MatSetSizes(pdipm->Jci_xb,pdipm->nci-pdipm->nh,pdipm->nx,PETSC_DECIDE,pdipm->Nx)); 1032 PetscCall(MatSetFromOptions(pdipm->Jci_xb)); 1033 PetscCall(MatSeqAIJSetPreallocation(pdipm->Jci_xb,1,NULL)); 1034 PetscCall(MatMPIAIJSetPreallocation(pdipm->Jci_xb,1,NULL,1,NULL)); 1035 1036 PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,&Jcrend)); 1037 offset = Jcrstart; 1038 if (pdipm->Nxub) { 1039 /* Add xub to Jci_xb */ 1040 PetscCall(ISGetIndices(pdipm->isxub,&cols)); 1041 k = 0; 1042 for (row = offset; row < offset + pdipm->nxub; row++) { 1043 PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES)); 1044 k++; 1045 } 1046 PetscCall(ISRestoreIndices(pdipm->isxub, &cols)); 1047 } 1048 1049 if (pdipm->Nxlb) { 1050 /* Add xlb to Jci_xb */ 1051 PetscCall(ISGetIndices(pdipm->isxlb,&cols)); 1052 k = 0; 1053 offset += pdipm->nxub; 1054 for (row = offset; row < offset + pdipm->nxlb; row++) { 1055 PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&one,INSERT_VALUES)); 1056 k++; 1057 } 1058 PetscCall(ISRestoreIndices(pdipm->isxlb, &cols)); 1059 } 1060 1061 /* Add xbox to Jci_xb */ 1062 if (pdipm->Nxbox) { 1063 PetscCall(ISGetIndices(pdipm->isxbox,&cols)); 1064 k = 0; 1065 offset += pdipm->nxlb; 1066 for (row = offset; row < offset + pdipm->nxbox; row++) { 1067 PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES)); 1068 tmp = row + pdipm->nxbox; 1069 PetscCall(MatSetValues(pdipm->Jci_xb,1,&tmp,1,cols+k,&one,INSERT_VALUES)); 1070 k++; 1071 } 1072 PetscCall(ISRestoreIndices(pdipm->isxbox, &cols)); 1073 } 1074 1075 PetscCall(MatAssemblyBegin(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY)); 1076 PetscCall(MatAssemblyEnd(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY)); 1077 /* PetscCall(MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD)); */ 1078 1079 /* (6) Set up ISs for PC Fieldsplit */ 1080 if (pdipm->solve_reduced_kkt) { 1081 PetscCall(PetscMalloc2(pdipm->nx+pdipm->nce,&xa,2*pdipm->nci,&xb)); 1082 for (i=0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i; 1083 for (i=0; i < 2*pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i; 1084 1085 PetscCall(ISCreateGeneral(comm,pdipm->nx+pdipm->nce,xa,PETSC_OWN_POINTER,&pdipm->is1)); 1086 PetscCall(ISCreateGeneral(comm,2*pdipm->nci,xb,PETSC_OWN_POINTER,&pdipm->is2)); 1087 } 1088 1089 /* (7) Gather offsets from all processes */ 1090 PetscCall(PetscMalloc1(size,&pdipm->nce_all)); 1091 1092 /* Get rstart of KKT matrix */ 1093 PetscCallMPI(MPI_Scan(&pdipm->n,&rstart,1,MPIU_INT,MPI_SUM,comm)); 1094 rstart -= pdipm->n; 1095 1096 PetscCallMPI(MPI_Allgather(&pdipm->nce,1,MPIU_INT,pdipm->nce_all,1,MPIU_INT,comm)); 1097 1098 PetscCall(PetscMalloc3(size,&ng_all,size,&nh_all,size,&Jranges)); 1099 PetscCallMPI(MPI_Allgather(&rstart,1,MPIU_INT,Jranges,1,MPIU_INT,comm)); 1100 PetscCallMPI(MPI_Allgather(&pdipm->nh,1,MPIU_INT,nh_all,1,MPIU_INT,comm)); 1101 PetscCallMPI(MPI_Allgather(&pdipm->ng,1,MPIU_INT,ng_all,1,MPIU_INT,comm)); 1102 1103 PetscCall(MatGetOwnershipRanges(tao->hessian,&rranges)); 1104 PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&cranges)); 1105 1106 if (pdipm->Ng) { 1107 PetscCall(TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre)); 1108 PetscCall(MatTranspose(tao->jacobian_equality,MAT_INITIAL_MATRIX,&pdipm->jac_equality_trans)); 1109 } 1110 if (pdipm->Nh) { 1111 PetscCall(TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre)); 1112 PetscCall(MatTranspose(tao->jacobian_inequality,MAT_INITIAL_MATRIX,&pdipm->jac_inequality_trans)); 1113 } 1114 1115 /* Count dnz,onz for preallocation of KKT matrix */ 1116 nce_all = pdipm->nce_all; 1117 1118 if (pdipm->Nxfixed) { 1119 PetscCall(MatTranspose(pdipm->Jce_xfixed,MAT_INITIAL_MATRIX,&Jce_xfixed_trans)); 1120 } 1121 PetscCall(MatTranspose(pdipm->Jci_xb,MAT_INITIAL_MATRIX,&Jci_xb_trans)); 1122 1123 MatPreallocateBegin(comm,pdipm->n,pdipm->n,dnz,onz); 1124 1125 /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */ 1126 PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao,pdipm->x)); 1127 PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 1128 1129 /* Insert tao->hessian */ 1130 PetscCall(MatGetOwnershipRange(tao->hessian,&rjstart,NULL)); 1131 for (i=0; i<pdipm->nx; i++) { 1132 row = rstart + i; 1133 1134 PetscCall(MatGetRow(tao->hessian,i+rjstart,&nc,&aj,NULL)); 1135 proc = 0; 1136 for (j=0; j < nc; j++) { 1137 while (aj[j] >= cranges[proc+1]) proc++; 1138 col = aj[j] - cranges[proc] + Jranges[proc]; 1139 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1140 } 1141 PetscCall(MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,NULL)); 1142 1143 if (pdipm->ng) { 1144 /* Insert grad g' */ 1145 PetscCall(MatGetRow(pdipm->jac_equality_trans,i+rjstart,&nc,&aj,NULL)); 1146 PetscCall(MatGetOwnershipRanges(tao->jacobian_equality,&ranges)); 1147 proc = 0; 1148 for (j=0; j < nc; j++) { 1149 /* find row ownership of */ 1150 while (aj[j] >= ranges[proc+1]) proc++; 1151 nx_all = rranges[proc+1] - rranges[proc]; 1152 col = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 1153 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1154 } 1155 PetscCall(MatRestoreRow(pdipm->jac_equality_trans,i+rjstart,&nc,&aj,NULL)); 1156 } 1157 1158 /* Insert Jce_xfixed^T' */ 1159 if (pdipm->nxfixed) { 1160 PetscCall(MatGetRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL)); 1161 PetscCall(MatGetOwnershipRanges(pdipm->Jce_xfixed,&ranges)); 1162 proc = 0; 1163 for (j=0; j < nc; j++) { 1164 /* find row ownership of */ 1165 while (aj[j] >= ranges[proc+1]) proc++; 1166 nx_all = rranges[proc+1] - rranges[proc]; 1167 col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc]; 1168 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1169 } 1170 PetscCall(MatRestoreRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL)); 1171 } 1172 1173 if (pdipm->nh) { 1174 /* Insert -grad h' */ 1175 PetscCall(MatGetRow(pdipm->jac_inequality_trans,i+rjstart,&nc,&aj,NULL)); 1176 PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality,&ranges)); 1177 proc = 0; 1178 for (j=0; j < nc; j++) { 1179 /* find row ownership of */ 1180 while (aj[j] >= ranges[proc+1]) proc++; 1181 nx_all = rranges[proc+1] - rranges[proc]; 1182 col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 1183 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1184 } 1185 PetscCall(MatRestoreRow(pdipm->jac_inequality_trans,i+rjstart,&nc,&aj,NULL)); 1186 } 1187 1188 /* Insert Jci_xb^T' */ 1189 PetscCall(MatGetRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL)); 1190 PetscCall(MatGetOwnershipRanges(pdipm->Jci_xb,&ranges)); 1191 proc = 0; 1192 for (j=0; j < nc; j++) { 1193 /* find row ownership of */ 1194 while (aj[j] >= ranges[proc+1]) proc++; 1195 nx_all = rranges[proc+1] - rranges[proc]; 1196 col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc]; 1197 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1198 } 1199 PetscCall(MatRestoreRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL)); 1200 } 1201 1202 /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */ 1203 if (pdipm->Ng) { 1204 PetscCall(MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL)); 1205 for (i=0; i < pdipm->ng; i++) { 1206 row = rstart + pdipm->off_lambdae + i; 1207 1208 PetscCall(MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL)); 1209 proc = 0; 1210 for (j=0; j < nc; j++) { 1211 while (aj[j] >= cranges[proc+1]) proc++; 1212 col = aj[j] - cranges[proc] + Jranges[proc]; 1213 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); /* grad g */ 1214 } 1215 PetscCall(MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL)); 1216 } 1217 } 1218 /* Jce_xfixed */ 1219 if (pdipm->Nxfixed) { 1220 PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL)); 1221 for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1222 row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1223 1224 PetscCall(MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL)); 1225 PetscCheck(nc == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1226 1227 proc = 0; 1228 j = 0; 1229 while (cols[j] >= cranges[proc+1]) proc++; 1230 col = cols[j] - cranges[proc] + Jranges[proc]; 1231 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1232 PetscCall(MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL)); 1233 } 1234 } 1235 1236 /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */ 1237 if (pdipm->Nh) { 1238 PetscCall(MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL)); 1239 for (i=0; i < pdipm->nh; i++) { 1240 row = rstart + pdipm->off_lambdai + i; 1241 1242 PetscCall(MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL)); 1243 proc = 0; 1244 for (j=0; j < nc; j++) { 1245 while (aj[j] >= cranges[proc+1]) proc++; 1246 col = aj[j] - cranges[proc] + Jranges[proc]; 1247 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); /* grad h */ 1248 } 1249 PetscCall(MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL)); 1250 } 1251 /* I */ 1252 for (i=0; i < pdipm->nh; i++) { 1253 row = rstart + pdipm->off_lambdai + i; 1254 col = rstart + pdipm->off_z + i; 1255 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1256 } 1257 } 1258 1259 /* Jci_xb */ 1260 PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL)); 1261 for (i=0; i < (pdipm->nci - pdipm->nh); i++) { 1262 row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1263 1264 PetscCall(MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL)); 1265 PetscCheck(nc == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1266 proc = 0; 1267 for (j=0; j < nc; j++) { 1268 while (cols[j] >= cranges[proc+1]) proc++; 1269 col = cols[j] - cranges[proc] + Jranges[proc]; 1270 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1271 } 1272 PetscCall(MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL)); 1273 /* I */ 1274 col = rstart + pdipm->off_z + pdipm->nh + i; 1275 PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1276 } 1277 1278 /* 4-th Row block of KKT matrix: Z and Ci */ 1279 for (i=0; i < pdipm->nci; i++) { 1280 row = rstart + pdipm->off_z + i; 1281 cols1[0] = rstart + pdipm->off_lambdai + i; 1282 cols1[1] = row; 1283 PetscCall(MatPreallocateSet(row,2,cols1,dnz,onz)); 1284 } 1285 1286 /* diagonal entry */ 1287 for (i=0; i<pdipm->n; i++) dnz[i]++; /* diagonal entry */ 1288 1289 /* Create KKT matrix */ 1290 PetscCall(MatCreate(comm,&J)); 1291 PetscCall(MatSetSizes(J,pdipm->n,pdipm->n,PETSC_DECIDE,PETSC_DECIDE)); 1292 PetscCall(MatSetFromOptions(J)); 1293 PetscCall(MatSeqAIJSetPreallocation(J,0,dnz)); 1294 PetscCall(MatMPIAIJSetPreallocation(J,0,dnz,0,onz)); 1295 MatPreallocateEnd(dnz,onz); 1296 pdipm->K = J; 1297 1298 /* (8) Insert constant entries to K */ 1299 /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */ 1300 PetscCall(MatGetOwnershipRange(J,&rstart,&rend)); 1301 for (i=rstart; i<rend; i++) { 1302 PetscCall(MatSetValue(J,i,i,0.0,INSERT_VALUES)); 1303 } 1304 /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */ 1305 if (pdipm->kkt_pd) { 1306 for (i=0; i<pdipm->nh; i++) { 1307 row = rstart + i; 1308 PetscCall(MatSetValue(J,row,row,pdipm->deltaw,INSERT_VALUES)); 1309 } 1310 } 1311 1312 /* Row block of K: [ grad Ce, 0, 0, 0] */ 1313 if (pdipm->Nxfixed) { 1314 PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL)); 1315 for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1316 row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1317 1318 PetscCall(MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa)); 1319 proc = 0; 1320 for (j=0; j < nc; j++) { 1321 while (cols[j] >= cranges[proc+1]) proc++; 1322 col = cols[j] - cranges[proc] + Jranges[proc]; 1323 PetscCall(MatSetValue(J,row,col,aa[j],INSERT_VALUES)); /* grad Ce */ 1324 PetscCall(MatSetValue(J,col,row,aa[j],INSERT_VALUES)); /* grad Ce' */ 1325 } 1326 PetscCall(MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa)); 1327 } 1328 } 1329 1330 /* Row block of K: [ -grad Ci, 0, 0, I] */ 1331 PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL)); 1332 for (i=0; i < pdipm->nci - pdipm->nh; i++) { 1333 row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1334 1335 PetscCall(MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa)); 1336 proc = 0; 1337 for (j=0; j < nc; j++) { 1338 while (cols[j] >= cranges[proc+1]) proc++; 1339 col = cols[j] - cranges[proc] + Jranges[proc]; 1340 PetscCall(MatSetValue(J,col,row,-aa[j],INSERT_VALUES)); 1341 PetscCall(MatSetValue(J,row,col,-aa[j],INSERT_VALUES)); 1342 } 1343 PetscCall(MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa)); 1344 1345 col = rstart + pdipm->off_z + pdipm->nh + i; 1346 PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 1347 } 1348 1349 for (i=0; i < pdipm->nh; i++) { 1350 row = rstart + pdipm->off_lambdai + i; 1351 col = rstart + pdipm->off_z + i; 1352 PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 1353 } 1354 1355 /* Row block of K: [ 0, 0, I, ...] */ 1356 for (i=0; i < pdipm->nci; i++) { 1357 row = rstart + pdipm->off_z + i; 1358 col = rstart + pdipm->off_lambdai + i; 1359 PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 1360 } 1361 1362 if (pdipm->Nxfixed) { 1363 PetscCall(MatDestroy(&Jce_xfixed_trans)); 1364 } 1365 PetscCall(MatDestroy(&Jci_xb_trans)); 1366 PetscCall(PetscFree3(ng_all,nh_all,Jranges)); 1367 1368 /* (9) Set up nonlinear solver SNES */ 1369 PetscCall(SNESSetFunction(pdipm->snes,NULL,TaoSNESFunction_PDIPM,(void*)tao)); 1370 PetscCall(SNESSetJacobian(pdipm->snes,J,J,TaoSNESJacobian_PDIPM,(void*)tao)); 1371 1372 if (pdipm->solve_reduced_kkt) { 1373 PC pc; 1374 PetscCall(KSPGetPC(tao->ksp,&pc)); 1375 PetscCall(PCSetType(pc,PCFIELDSPLIT)); 1376 PetscCall(PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR)); 1377 PetscCall(PCFieldSplitSetIS(pc,"2",pdipm->is2)); 1378 PetscCall(PCFieldSplitSetIS(pc,"1",pdipm->is1)); 1379 } 1380 PetscCall(SNESSetFromOptions(pdipm->snes)); 1381 1382 /* (10) Setup PCPreSolve() for pdipm->solve_symmetric_kkt */ 1383 if (pdipm->solve_symmetric_kkt) { 1384 KSP ksp; 1385 PC pc; 1386 PetscBool isCHOL; 1387 PetscCall(SNESGetKSP(pdipm->snes,&ksp)); 1388 PetscCall(KSPGetPC(ksp,&pc)); 1389 PetscCall(PCSetPreSolve(pc,PCPreSolve_PDIPM)); 1390 1391 PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL)); 1392 if (isCHOL) { 1393 Mat Factor; 1394 PetscBool isMUMPS; 1395 PetscCall(PCFactorGetMatrix(pc,&Factor)); 1396 PetscCall(PetscObjectTypeCompare((PetscObject)Factor,"mumps",&isMUMPS)); 1397 if (isMUMPS) { /* must set mumps ICNTL(13)=1 and ICNTL(24)=1 to call MatGetInertia() */ 1398 #if defined(PETSC_HAVE_MUMPS) 1399 PetscCall(MatMumpsSetIcntl(Factor,24,1)); /* detection of null pivot rows */ 1400 if (size > 1) { 1401 PetscCall(MatMumpsSetIcntl(Factor,13,1)); /* parallelism of the root node (enable ScaLAPACK) and its splitting */ 1402 } 1403 #else 1404 SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Requires external package MUMPS"); 1405 #endif 1406 } 1407 } 1408 } 1409 PetscFunctionReturn(0); 1410 } 1411 1412 /* 1413 TaoDestroy_PDIPM - Destroys the pdipm object 1414 1415 Input: 1416 full pdipm 1417 1418 Output: 1419 Destroyed pdipm 1420 */ 1421 PetscErrorCode TaoDestroy_PDIPM(Tao tao) 1422 { 1423 TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1424 1425 PetscFunctionBegin; 1426 /* Freeing Vectors assocaiated with KKT (X) */ 1427 PetscCall(VecDestroy(&pdipm->x)); /* Solution x */ 1428 PetscCall(VecDestroy(&pdipm->lambdae)); /* Equality constraints lagrangian multiplier*/ 1429 PetscCall(VecDestroy(&pdipm->lambdai)); /* Inequality constraints lagrangian multiplier*/ 1430 PetscCall(VecDestroy(&pdipm->z)); /* Slack variables */ 1431 PetscCall(VecDestroy(&pdipm->X)); /* Big KKT system vector [x; lambdae; lambdai; z] */ 1432 1433 /* work vectors */ 1434 PetscCall(VecDestroy(&pdipm->lambdae_xfixed)); 1435 PetscCall(VecDestroy(&pdipm->lambdai_xb)); 1436 1437 /* Legrangian equality and inequality Vec */ 1438 PetscCall(VecDestroy(&pdipm->ce)); /* Vec of equality constraints */ 1439 PetscCall(VecDestroy(&pdipm->ci)); /* Vec of inequality constraints */ 1440 1441 /* Matrices */ 1442 PetscCall(MatDestroy(&pdipm->Jce_xfixed)); 1443 PetscCall(MatDestroy(&pdipm->Jci_xb)); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */ 1444 PetscCall(MatDestroy(&pdipm->K)); 1445 1446 /* Index Sets */ 1447 if (pdipm->Nxub) { 1448 PetscCall(ISDestroy(&pdipm->isxub)); /* Finite upper bound only -inf < x < ub */ 1449 } 1450 1451 if (pdipm->Nxlb) { 1452 PetscCall(ISDestroy(&pdipm->isxlb)); /* Finite lower bound only lb <= x < inf */ 1453 } 1454 1455 if (pdipm->Nxfixed) { 1456 PetscCall(ISDestroy(&pdipm->isxfixed)); /* Fixed variables lb = x = ub */ 1457 } 1458 1459 if (pdipm->Nxbox) { 1460 PetscCall(ISDestroy(&pdipm->isxbox)); /* Boxed variables lb <= x <= ub */ 1461 } 1462 1463 if (pdipm->Nxfree) { 1464 PetscCall(ISDestroy(&pdipm->isxfree)); /* Free variables -inf <= x <= inf */ 1465 } 1466 1467 if (pdipm->solve_reduced_kkt) { 1468 PetscCall(ISDestroy(&pdipm->is1)); 1469 PetscCall(ISDestroy(&pdipm->is2)); 1470 } 1471 1472 /* SNES */ 1473 PetscCall(SNESDestroy(&pdipm->snes)); /* Nonlinear solver */ 1474 PetscCall(PetscFree(pdipm->nce_all)); 1475 PetscCall(MatDestroy(&pdipm->jac_equality_trans)); 1476 PetscCall(MatDestroy(&pdipm->jac_inequality_trans)); 1477 1478 /* Destroy pdipm */ 1479 PetscCall(PetscFree(tao->data)); /* Holding locations of pdipm */ 1480 1481 /* Destroy Dual */ 1482 PetscCall(VecDestroy(&tao->DE)); /* equality dual */ 1483 PetscCall(VecDestroy(&tao->DI)); /* dinequality dual */ 1484 PetscFunctionReturn(0); 1485 } 1486 1487 PetscErrorCode TaoSetFromOptions_PDIPM(PetscOptionItems *PetscOptionsObject,Tao tao) 1488 { 1489 TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1490 1491 PetscFunctionBegin; 1492 PetscOptionsHeadBegin(PetscOptionsObject,"PDIPM method for constrained optimization"); 1493 PetscCall(PetscOptionsReal("-tao_pdipm_push_init_slack","parameter to push initial slack variables away from bounds",NULL,pdipm->push_init_slack,&pdipm->push_init_slack,NULL)); 1494 PetscCall(PetscOptionsReal("-tao_pdipm_push_init_lambdai","parameter to push initial (inequality) dual variables away from bounds",NULL,pdipm->push_init_lambdai,&pdipm->push_init_lambdai,NULL)); 1495 PetscCall(PetscOptionsBool("-tao_pdipm_solve_reduced_kkt","Solve reduced KKT system using Schur-complement",NULL,pdipm->solve_reduced_kkt,&pdipm->solve_reduced_kkt,NULL)); 1496 PetscCall(PetscOptionsReal("-tao_pdipm_mu_update_factor","Update scalar for barrier parameter (mu) update",NULL,pdipm->mu_update_factor,&pdipm->mu_update_factor,NULL)); 1497 PetscCall(PetscOptionsBool("-tao_pdipm_symmetric_kkt","Solve non reduced symmetric KKT system",NULL,pdipm->solve_symmetric_kkt,&pdipm->solve_symmetric_kkt,NULL)); 1498 PetscCall(PetscOptionsBool("-tao_pdipm_kkt_shift_pd","Add shifts to make KKT matrix positive definite",NULL,pdipm->kkt_pd,&pdipm->kkt_pd,NULL)); 1499 PetscOptionsHeadEnd(); 1500 PetscFunctionReturn(0); 1501 } 1502 1503 /*MC 1504 TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization. 1505 1506 Option Database Keys: 1507 + -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0) 1508 . -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0) 1509 . -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0) 1510 . -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system 1511 - -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite 1512 1513 Level: beginner 1514 M*/ 1515 PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao) 1516 { 1517 TAO_PDIPM *pdipm; 1518 1519 PetscFunctionBegin; 1520 tao->ops->setup = TaoSetup_PDIPM; 1521 tao->ops->solve = TaoSolve_PDIPM; 1522 tao->ops->setfromoptions = TaoSetFromOptions_PDIPM; 1523 tao->ops->view = TaoView_PDIPM; 1524 tao->ops->destroy = TaoDestroy_PDIPM; 1525 1526 PetscCall(PetscNewLog(tao,&pdipm)); 1527 tao->data = (void*)pdipm; 1528 1529 pdipm->nx = pdipm->Nx = 0; 1530 pdipm->nxfixed = pdipm->Nxfixed = 0; 1531 pdipm->nxlb = pdipm->Nxlb = 0; 1532 pdipm->nxub = pdipm->Nxub = 0; 1533 pdipm->nxbox = pdipm->Nxbox = 0; 1534 pdipm->nxfree = pdipm->Nxfree = 0; 1535 1536 pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0; 1537 pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0; 1538 pdipm->n = pdipm->N = 0; 1539 pdipm->mu = 1.0; 1540 pdipm->mu_update_factor = 0.1; 1541 1542 pdipm->deltaw = 0.0; 1543 pdipm->lastdeltaw = 3*1.e-4; 1544 pdipm->deltac = 0.0; 1545 pdipm->kkt_pd = PETSC_FALSE; 1546 1547 pdipm->push_init_slack = 1.0; 1548 pdipm->push_init_lambdai = 1.0; 1549 pdipm->solve_reduced_kkt = PETSC_FALSE; 1550 pdipm->solve_symmetric_kkt = PETSC_TRUE; 1551 1552 /* Override default settings (unless already changed) */ 1553 if (!tao->max_it_changed) tao->max_it = 200; 1554 if (!tao->max_funcs_changed) tao->max_funcs = 500; 1555 1556 PetscCall(SNESCreate(((PetscObject)tao)->comm,&pdipm->snes)); 1557 PetscCall(SNESSetOptionsPrefix(pdipm->snes,tao->hdr.prefix)); 1558 PetscCall(SNESGetKSP(pdipm->snes,&tao->ksp)); 1559 PetscCall(PetscObjectReference((PetscObject)tao->ksp)); 1560 PetscCall(KSPSetApplicationContext(tao->ksp,(void *)tao)); 1561 PetscFunctionReturn(0); 1562 } 1563