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