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