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: `TAOPDIPM`, `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: `TAOPDIPM`, `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: `TAOPDIPM`, `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(MPIU_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 `TAOPDIPM` 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 `VecPlaceArray()`/`VecResetArray()` 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 - matrix to build the preconditioner from 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 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 PCPostSetup_PDIPM -- called when the KKT matrix is Cholesky factored for the preconditioner. Checks the inertia of Cholesky factor of the 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 PCPostSetUp_PDIPM(PC pc) 625 { 626 Tao tao; 627 TAO_PDIPM *pdipm; 628 Vec X; 629 SNES snes; 630 KSP ksp; 631 Mat Factor; 632 PetscBool isCHOL; 633 PetscInt nneg, nzero, npos; 634 635 PetscFunctionBegin; 636 PetscCall(PCGetApplicationContext(pc, &tao)); 637 pdipm = (TAO_PDIPM *)tao->data; 638 X = pdipm->X; 639 snes = pdipm->snes; 640 641 /* Get the inertia of Cholesky factor */ 642 PetscCall(SNESGetKSP(snes, &ksp)); 643 PetscCall(KSPGetPC(ksp, &pc)); 644 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCCHOLESKY, &isCHOL)); 645 if (!isCHOL) PetscFunctionReturn(PETSC_SUCCESS); 646 647 PetscCall(PCFactorGetMatrix(pc, &Factor)); 648 PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos)); 649 650 if (npos < pdipm->Nx + pdipm->Nci) { 651 pdipm->deltaw = PetscMax(pdipm->lastdeltaw / 3, 1.e-4 * PETSC_MACHINE_EPSILON); 652 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)); 653 PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao)); 654 PetscCall(PCSetPostSetUp(pc, NULL)); 655 PetscCall(PCSetUp(pc)); 656 PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos)); 657 658 if (npos < pdipm->Nx + pdipm->Nci) { 659 pdipm->deltaw = pdipm->lastdeltaw; /* in case reduction update does not help, this prevents that step from impacting increasing update */ 660 while (npos < pdipm->Nx + pdipm->Nci && pdipm->deltaw <= 1. / PETSC_SMALL) { /* increase deltaw */ 661 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)); 662 pdipm->deltaw = PetscMin(8 * pdipm->deltaw, PetscPowReal(10, 20)); 663 PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao)); 664 PetscCall(PCSetUp(pc)); 665 PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos)); 666 } 667 668 PetscCheck(pdipm->deltaw < 1. / PETSC_SMALL, PetscObjectComm((PetscObject)tao), PETSC_ERR_CONV_FAILED, "Reached maximum delta w will not converge, try different initial x0"); 669 670 PetscCall(PetscInfo(tao, "Updated deltaw %g\n", (double)pdipm->deltaw)); 671 pdipm->lastdeltaw = pdipm->deltaw; 672 pdipm->deltaw = 0.0; 673 } 674 } 675 676 if (nzero) { /* Jacobian is singular */ 677 if (pdipm->deltac == 0.0) { 678 pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON; 679 } else { 680 pdipm->deltac = pdipm->deltac * PetscPowReal(pdipm->mu, .25); 681 } 682 PetscCall(PetscInfo(tao, "Updated deltac=%g, MatInertia: nneg %" PetscInt_FMT ", nzero %" PetscInt_FMT "(!=0), npos %" PetscInt_FMT "\n", (double)pdipm->deltac, nneg, nzero, npos)); 683 PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao)); 684 PetscCall(PCSetPostSetUp(pc, NULL)); 685 PetscCall(PCSetUp(pc)); 686 PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos)); 687 } 688 PetscCall(PCSetPostSetUp(pc, PCPostSetUp_PDIPM)); 689 PetscFunctionReturn(PETSC_SUCCESS); 690 } 691 692 /* 693 SNESLineSearch_PDIPM - Custom line search used with PDIPM. 694 695 Collective 696 697 Notes: 698 This routine employs a simple backtracking line-search to keep 699 the slack variables (z) and inequality constraints Lagrange multipliers 700 (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars 701 alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the 702 slack variables are updated as X = X - alpha_d*dx. The constraint multipliers 703 are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu 704 is also updated as mu = mu + z'lambdai/Nci 705 */ 706 static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch, void *ctx) 707 { 708 Tao tao = (Tao)ctx; 709 TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 710 SNES snes; 711 Vec X, F, Y; 712 PetscInt i, iter; 713 PetscReal alpha_p = 1.0, alpha_d = 1.0, alpha[4]; 714 PetscScalar *Xarr, *z, *lambdai, dot, *taosolarr; 715 const PetscScalar *dXarr, *dz, *dlambdai; 716 717 PetscFunctionBegin; 718 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 719 PetscCall(SNESGetIterationNumber(snes, &iter)); 720 721 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 722 PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, NULL, NULL)); 723 724 PetscCall(VecGetArrayWrite(X, &Xarr)); 725 PetscCall(VecGetArrayRead(Y, &dXarr)); 726 z = Xarr + pdipm->off_z; 727 dz = dXarr + pdipm->off_z; 728 for (i = 0; i < pdipm->nci; i++) { 729 if (z[i] - dz[i] < 0.0) alpha_p = PetscMin(alpha_p, 0.9999 * z[i] / dz[i]); 730 } 731 732 lambdai = Xarr + pdipm->off_lambdai; 733 dlambdai = dXarr + pdipm->off_lambdai; 734 735 for (i = 0; i < pdipm->nci; i++) { 736 if (lambdai[i] - dlambdai[i] < 0.0) alpha_d = PetscMin(0.9999 * lambdai[i] / dlambdai[i], alpha_d); 737 } 738 739 alpha[0] = alpha_p; 740 alpha[1] = alpha_d; 741 PetscCall(VecRestoreArrayRead(Y, &dXarr)); 742 PetscCall(VecRestoreArrayWrite(X, &Xarr)); 743 744 /* alpha = min(alpha) over all processes */ 745 PetscCallMPI(MPIU_Allreduce(alpha, alpha + 2, 2, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)tao))); 746 747 alpha_p = alpha[2]; 748 alpha_d = alpha[3]; 749 750 /* X = X - alpha * Y */ 751 PetscCall(VecGetArrayWrite(X, &Xarr)); 752 PetscCall(VecGetArrayRead(Y, &dXarr)); 753 for (i = 0; i < pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i]; 754 for (i = 0; i < pdipm->nce; i++) Xarr[i + pdipm->off_lambdae] -= alpha_d * dXarr[i + pdipm->off_lambdae]; 755 756 for (i = 0; i < pdipm->nci; i++) { 757 Xarr[i + pdipm->off_lambdai] -= alpha_d * dXarr[i + pdipm->off_lambdai]; 758 Xarr[i + pdipm->off_z] -= alpha_p * dXarr[i + pdipm->off_z]; 759 } 760 PetscCall(VecGetArrayWrite(tao->solution, &taosolarr)); 761 PetscCall(PetscMemcpy(taosolarr, Xarr, pdipm->nx * sizeof(PetscScalar))); 762 PetscCall(VecRestoreArrayWrite(tao->solution, &taosolarr)); 763 764 PetscCall(VecRestoreArrayWrite(X, &Xarr)); 765 PetscCall(VecRestoreArrayRead(Y, &dXarr)); 766 767 /* Update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */ 768 if (pdipm->z) PetscCall(VecDot(pdipm->z, pdipm->lambdai, &dot)); 769 else dot = 0.0; 770 771 /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu) */ 772 pdipm->mu = pdipm->mu_update_factor * dot / pdipm->Nci; 773 774 /* Update F; get tao->residual and tao->cnorm */ 775 PetscCall(TaoSNESFunction_PDIPM_residual(snes, X, F, (void *)tao)); 776 777 tao->niter++; 778 PetscCall(TaoLogConvergenceHistory(tao, pdipm->obj, tao->residual, tao->cnorm, tao->niter)); 779 PetscCall(TaoMonitor(tao, tao->niter, pdipm->obj, tao->residual, tao->cnorm, pdipm->mu)); 780 781 PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 782 if (tao->reason) PetscCall(SNESSetConvergedReason(snes, SNES_CONVERGED_FNORM_ABS)); 783 PetscFunctionReturn(PETSC_SUCCESS); 784 } 785 786 static PetscErrorCode TaoSolve_PDIPM(Tao tao) 787 { 788 TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 789 SNESLineSearch linesearch; /* SNESLineSearch context */ 790 Vec dummy; 791 792 PetscFunctionBegin; 793 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"); 794 795 /* Initialize all variables */ 796 PetscCall(TaoPDIPMInitializeSolution(tao)); 797 798 /* Set linesearch */ 799 PetscCall(SNESGetLineSearch(pdipm->snes, &linesearch)); 800 PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSHELL)); 801 PetscCall(SNESLineSearchShellSetApply(linesearch, SNESLineSearch_PDIPM, tao)); 802 PetscCall(SNESLineSearchSetFromOptions(linesearch)); 803 804 tao->reason = TAO_CONTINUE_ITERATING; 805 806 /* -tao_monitor for iteration 0 and check convergence */ 807 PetscCall(VecDuplicate(pdipm->X, &dummy)); 808 PetscCall(TaoSNESFunction_PDIPM_residual(pdipm->snes, pdipm->X, dummy, (void *)tao)); 809 810 PetscCall(TaoLogConvergenceHistory(tao, pdipm->obj, tao->residual, tao->cnorm, tao->niter)); 811 PetscCall(TaoMonitor(tao, tao->niter, pdipm->obj, tao->residual, tao->cnorm, pdipm->mu)); 812 PetscCall(VecDestroy(&dummy)); 813 PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 814 if (tao->reason) PetscCall(SNESSetConvergedReason(pdipm->snes, SNES_CONVERGED_FNORM_ABS)); 815 816 while (tao->reason == TAO_CONTINUE_ITERATING) { 817 SNESConvergedReason reason; 818 PetscCall(SNESSolve(pdipm->snes, NULL, pdipm->X)); 819 820 /* Check SNES convergence */ 821 PetscCall(SNESGetConvergedReason(pdipm->snes, &reason)); 822 if (reason < 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes), "SNES solve did not converged due to reason %s\n", SNESConvergedReasons[reason])); 823 824 /* Check TAO convergence */ 825 PetscCheck(!PetscIsInfOrNanReal(pdipm->obj), PETSC_COMM_SELF, PETSC_ERR_SUP, "User-provided compute function generated Inf or NaN"); 826 } 827 PetscFunctionReturn(PETSC_SUCCESS); 828 } 829 830 static PetscErrorCode TaoView_PDIPM(Tao tao, PetscViewer viewer) 831 { 832 TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 833 834 PetscFunctionBegin; 835 tao->constrained = PETSC_TRUE; 836 PetscCall(PetscViewerASCIIPushTab(viewer)); 837 PetscCall(PetscViewerASCIIPrintf(viewer, "Number of prime=%" PetscInt_FMT ", Number of dual=%" PetscInt_FMT "\n", pdipm->Nx + pdipm->Nci, pdipm->Nce + pdipm->Nci)); 838 if (pdipm->kkt_pd) PetscCall(PetscViewerASCIIPrintf(viewer, "KKT shifts deltaw=%g, deltac=%g\n", (double)pdipm->deltaw, (double)pdipm->deltac)); 839 PetscCall(PetscViewerASCIIPopTab(viewer)); 840 PetscFunctionReturn(PETSC_SUCCESS); 841 } 842 843 static PetscErrorCode TaoSetup_PDIPM(Tao tao) 844 { 845 TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 846 MPI_Comm comm; 847 PetscMPIInt size; 848 PetscInt row, col, Jcrstart, Jcrend, k, tmp, nc, proc, *nh_all, *ng_all; 849 PetscInt offset, *xa, *xb, i, j, rstart, rend; 850 PetscScalar one = 1.0, neg_one = -1.0; 851 const PetscInt *cols, *rranges, *cranges, *aj, *ranges; 852 const PetscScalar *aa, *Xarr; 853 Mat J; 854 Mat Jce_xfixed_trans, Jci_xb_trans; 855 PetscInt *dnz, *onz, rjstart, nx_all, *nce_all, *Jranges, cols1[2]; 856 857 PetscFunctionBegin; 858 PetscCall(PetscObjectGetComm((PetscObject)tao, &comm)); 859 PetscCallMPI(MPI_Comm_size(comm, &size)); 860 861 /* (1) Setup Bounds and create Tao vectors */ 862 PetscCall(TaoPDIPMSetUpBounds(tao)); 863 864 if (!tao->gradient) { 865 PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 866 PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 867 } 868 869 /* (2) Get sizes */ 870 /* Size of vector x - This is set by TaoSetSolution */ 871 PetscCall(VecGetSize(tao->solution, &pdipm->Nx)); 872 PetscCall(VecGetLocalSize(tao->solution, &pdipm->nx)); 873 874 /* Size of equality constraints and vectors */ 875 if (tao->constraints_equality) { 876 PetscCall(VecGetSize(tao->constraints_equality, &pdipm->Ng)); 877 PetscCall(VecGetLocalSize(tao->constraints_equality, &pdipm->ng)); 878 } else { 879 pdipm->ng = pdipm->Ng = 0; 880 } 881 882 pdipm->nce = pdipm->ng + pdipm->nxfixed; 883 pdipm->Nce = pdipm->Ng + pdipm->Nxfixed; 884 885 /* Size of inequality constraints and vectors */ 886 if (tao->constraints_inequality) { 887 PetscCall(VecGetSize(tao->constraints_inequality, &pdipm->Nh)); 888 PetscCall(VecGetLocalSize(tao->constraints_inequality, &pdipm->nh)); 889 } else { 890 pdipm->nh = pdipm->Nh = 0; 891 } 892 893 pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2 * pdipm->nxbox; 894 pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2 * pdipm->Nxbox; 895 896 /* Full size of the KKT system to be solved */ 897 pdipm->n = pdipm->nx + pdipm->nce + 2 * pdipm->nci; 898 pdipm->N = pdipm->Nx + pdipm->Nce + 2 * pdipm->Nci; 899 900 /* (3) Offsets for subvectors */ 901 pdipm->off_lambdae = pdipm->nx; 902 pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce; 903 pdipm->off_z = pdipm->off_lambdai + pdipm->nci; 904 905 /* (4) Create vectors and subvectors */ 906 /* Ce and Ci vectors */ 907 PetscCall(VecCreate(comm, &pdipm->ce)); 908 PetscCall(VecSetSizes(pdipm->ce, pdipm->nce, pdipm->Nce)); 909 PetscCall(VecSetFromOptions(pdipm->ce)); 910 911 PetscCall(VecCreate(comm, &pdipm->ci)); 912 PetscCall(VecSetSizes(pdipm->ci, pdipm->nci, pdipm->Nci)); 913 PetscCall(VecSetFromOptions(pdipm->ci)); 914 915 /* X=[x; lambdae; lambdai; z] for the big KKT system */ 916 PetscCall(VecCreate(comm, &pdipm->X)); 917 PetscCall(VecSetSizes(pdipm->X, pdipm->n, pdipm->N)); 918 PetscCall(VecSetFromOptions(pdipm->X)); 919 920 /* Subvectors; they share local arrays with X */ 921 PetscCall(VecGetArrayRead(pdipm->X, &Xarr)); 922 /* x shares local array with X.x */ 923 if (pdipm->Nx) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nx, pdipm->Nx, Xarr, &pdipm->x)); 924 925 /* lambdae shares local array with X.lambdae */ 926 if (pdipm->Nce) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nce, pdipm->Nce, Xarr + pdipm->off_lambdae, &pdipm->lambdae)); 927 928 /* tao->DE shares local array with X.lambdae_g */ 929 if (pdipm->Ng) { 930 PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->ng, pdipm->Ng, Xarr + pdipm->off_lambdae, &tao->DE)); 931 932 PetscCall(VecCreate(comm, &pdipm->lambdae_xfixed)); 933 PetscCall(VecSetSizes(pdipm->lambdae_xfixed, pdipm->nxfixed, PETSC_DECIDE)); 934 PetscCall(VecSetFromOptions(pdipm->lambdae_xfixed)); 935 } 936 937 if (pdipm->Nci) { 938 /* lambdai shares local array with X.lambdai */ 939 PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nci, pdipm->Nci, Xarr + pdipm->off_lambdai, &pdipm->lambdai)); 940 941 /* z for slack variables; it shares local array with X.z */ 942 PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nci, pdipm->Nci, Xarr + pdipm->off_z, &pdipm->z)); 943 } 944 945 /* tao->DI which shares local array with X.lambdai_h */ 946 if (pdipm->Nh) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nh, pdipm->Nh, Xarr + pdipm->off_lambdai, &tao->DI)); 947 PetscCall(VecCreate(comm, &pdipm->lambdai_xb)); 948 PetscCall(VecSetSizes(pdipm->lambdai_xb, pdipm->nci - pdipm->nh, PETSC_DECIDE)); 949 PetscCall(VecSetFromOptions(pdipm->lambdai_xb)); 950 951 PetscCall(VecRestoreArrayRead(pdipm->X, &Xarr)); 952 953 /* (5) Create Jacobians Jce_xfixed and Jci */ 954 /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */ 955 if (pdipm->Nxfixed) { 956 /* Create Jce_xfixed */ 957 PetscCall(MatCreate(comm, &pdipm->Jce_xfixed)); 958 PetscCall(MatSetSizes(pdipm->Jce_xfixed, pdipm->nxfixed, pdipm->nx, PETSC_DECIDE, pdipm->Nx)); 959 PetscCall(MatSetFromOptions(pdipm->Jce_xfixed)); 960 PetscCall(MatSeqAIJSetPreallocation(pdipm->Jce_xfixed, 1, NULL)); 961 PetscCall(MatMPIAIJSetPreallocation(pdipm->Jce_xfixed, 1, NULL, 1, NULL)); 962 963 PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, &Jcrend)); 964 PetscCall(ISGetIndices(pdipm->isxfixed, &cols)); 965 k = 0; 966 for (row = Jcrstart; row < Jcrend; row++) { 967 PetscCall(MatSetValues(pdipm->Jce_xfixed, 1, &row, 1, cols + k, &one, INSERT_VALUES)); 968 k++; 969 } 970 PetscCall(ISRestoreIndices(pdipm->isxfixed, &cols)); 971 PetscCall(MatAssemblyBegin(pdipm->Jce_xfixed, MAT_FINAL_ASSEMBLY)); 972 PetscCall(MatAssemblyEnd(pdipm->Jce_xfixed, MAT_FINAL_ASSEMBLY)); 973 } 974 975 /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */ 976 PetscCall(MatCreate(comm, &pdipm->Jci_xb)); 977 PetscCall(MatSetSizes(pdipm->Jci_xb, pdipm->nci - pdipm->nh, pdipm->nx, PETSC_DECIDE, pdipm->Nx)); 978 PetscCall(MatSetFromOptions(pdipm->Jci_xb)); 979 PetscCall(MatSeqAIJSetPreallocation(pdipm->Jci_xb, 1, NULL)); 980 PetscCall(MatMPIAIJSetPreallocation(pdipm->Jci_xb, 1, NULL, 1, NULL)); 981 982 PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, &Jcrend)); 983 offset = Jcrstart; 984 if (pdipm->Nxub) { 985 /* Add xub to Jci_xb */ 986 PetscCall(ISGetIndices(pdipm->isxub, &cols)); 987 k = 0; 988 for (row = offset; row < offset + pdipm->nxub; row++) { 989 PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &neg_one, INSERT_VALUES)); 990 k++; 991 } 992 PetscCall(ISRestoreIndices(pdipm->isxub, &cols)); 993 } 994 995 if (pdipm->Nxlb) { 996 /* Add xlb to Jci_xb */ 997 PetscCall(ISGetIndices(pdipm->isxlb, &cols)); 998 k = 0; 999 offset += pdipm->nxub; 1000 for (row = offset; row < offset + pdipm->nxlb; row++) { 1001 PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &one, INSERT_VALUES)); 1002 k++; 1003 } 1004 PetscCall(ISRestoreIndices(pdipm->isxlb, &cols)); 1005 } 1006 1007 /* Add xbox to Jci_xb */ 1008 if (pdipm->Nxbox) { 1009 PetscCall(ISGetIndices(pdipm->isxbox, &cols)); 1010 k = 0; 1011 offset += pdipm->nxlb; 1012 for (row = offset; row < offset + pdipm->nxbox; row++) { 1013 PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &neg_one, INSERT_VALUES)); 1014 tmp = row + pdipm->nxbox; 1015 PetscCall(MatSetValues(pdipm->Jci_xb, 1, &tmp, 1, cols + k, &one, INSERT_VALUES)); 1016 k++; 1017 } 1018 PetscCall(ISRestoreIndices(pdipm->isxbox, &cols)); 1019 } 1020 1021 PetscCall(MatAssemblyBegin(pdipm->Jci_xb, MAT_FINAL_ASSEMBLY)); 1022 PetscCall(MatAssemblyEnd(pdipm->Jci_xb, MAT_FINAL_ASSEMBLY)); 1023 /* PetscCall(MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD)); */ 1024 1025 /* (6) Set up ISs for PC Fieldsplit */ 1026 if (pdipm->solve_reduced_kkt) { 1027 PetscCall(PetscMalloc2(pdipm->nx + pdipm->nce, &xa, 2 * pdipm->nci, &xb)); 1028 for (i = 0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i; 1029 for (i = 0; i < 2 * pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i; 1030 1031 PetscCall(ISCreateGeneral(comm, pdipm->nx + pdipm->nce, xa, PETSC_OWN_POINTER, &pdipm->is1)); 1032 PetscCall(ISCreateGeneral(comm, 2 * pdipm->nci, xb, PETSC_OWN_POINTER, &pdipm->is2)); 1033 } 1034 1035 /* (7) Gather offsets from all processes */ 1036 PetscCall(PetscMalloc1(size, &pdipm->nce_all)); 1037 1038 /* Get rstart of KKT matrix */ 1039 PetscCallMPI(MPI_Scan(&pdipm->n, &rstart, 1, MPIU_INT, MPI_SUM, comm)); 1040 rstart -= pdipm->n; 1041 1042 PetscCallMPI(MPI_Allgather(&pdipm->nce, 1, MPIU_INT, pdipm->nce_all, 1, MPIU_INT, comm)); 1043 1044 PetscCall(PetscMalloc3(size, &ng_all, size, &nh_all, size, &Jranges)); 1045 PetscCallMPI(MPI_Allgather(&rstart, 1, MPIU_INT, Jranges, 1, MPIU_INT, comm)); 1046 PetscCallMPI(MPI_Allgather(&pdipm->nh, 1, MPIU_INT, nh_all, 1, MPIU_INT, comm)); 1047 PetscCallMPI(MPI_Allgather(&pdipm->ng, 1, MPIU_INT, ng_all, 1, MPIU_INT, comm)); 1048 1049 PetscCall(MatGetOwnershipRanges(tao->hessian, &rranges)); 1050 PetscCall(MatGetOwnershipRangesColumn(tao->hessian, &cranges)); 1051 1052 if (pdipm->Ng) { 1053 PetscCall(TaoComputeJacobianEquality(tao, tao->solution, tao->jacobian_equality, tao->jacobian_equality_pre)); 1054 PetscCall(MatTranspose(tao->jacobian_equality, MAT_INITIAL_MATRIX, &pdipm->jac_equality_trans)); 1055 } 1056 if (pdipm->Nh) { 1057 PetscCall(TaoComputeJacobianInequality(tao, tao->solution, tao->jacobian_inequality, tao->jacobian_inequality_pre)); 1058 PetscCall(MatTranspose(tao->jacobian_inequality, MAT_INITIAL_MATRIX, &pdipm->jac_inequality_trans)); 1059 } 1060 1061 /* Count dnz,onz for preallocation of KKT matrix */ 1062 nce_all = pdipm->nce_all; 1063 1064 if (pdipm->Nxfixed) PetscCall(MatTranspose(pdipm->Jce_xfixed, MAT_INITIAL_MATRIX, &Jce_xfixed_trans)); 1065 PetscCall(MatTranspose(pdipm->Jci_xb, MAT_INITIAL_MATRIX, &Jci_xb_trans)); 1066 1067 MatPreallocateBegin(comm, pdipm->n, pdipm->n, dnz, onz); 1068 1069 /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */ 1070 PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao, pdipm->x)); 1071 PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre)); 1072 1073 /* Insert tao->hessian */ 1074 PetscCall(MatGetOwnershipRange(tao->hessian, &rjstart, NULL)); 1075 for (i = 0; i < pdipm->nx; i++) { 1076 row = rstart + i; 1077 1078 PetscCall(MatGetRow(tao->hessian, i + rjstart, &nc, &aj, NULL)); 1079 proc = 0; 1080 for (j = 0; j < nc; j++) { 1081 while (aj[j] >= cranges[proc + 1]) proc++; 1082 col = aj[j] - cranges[proc] + Jranges[proc]; 1083 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1084 } 1085 PetscCall(MatRestoreRow(tao->hessian, i + rjstart, &nc, &aj, NULL)); 1086 1087 if (pdipm->ng) { 1088 /* Insert grad g' */ 1089 PetscCall(MatGetRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, NULL)); 1090 PetscCall(MatGetOwnershipRanges(tao->jacobian_equality, &ranges)); 1091 proc = 0; 1092 for (j = 0; j < nc; j++) { 1093 /* find row ownership of */ 1094 while (aj[j] >= ranges[proc + 1]) proc++; 1095 nx_all = rranges[proc + 1] - rranges[proc]; 1096 col = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 1097 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1098 } 1099 PetscCall(MatRestoreRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, NULL)); 1100 } 1101 1102 /* Insert Jce_xfixed^T' */ 1103 if (pdipm->nxfixed) { 1104 PetscCall(MatGetRow(Jce_xfixed_trans, i + rjstart, &nc, &aj, NULL)); 1105 PetscCall(MatGetOwnershipRanges(pdipm->Jce_xfixed, &ranges)); 1106 proc = 0; 1107 for (j = 0; j < nc; j++) { 1108 /* find row ownership of */ 1109 while (aj[j] >= ranges[proc + 1]) proc++; 1110 nx_all = rranges[proc + 1] - rranges[proc]; 1111 col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc]; 1112 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1113 } 1114 PetscCall(MatRestoreRow(Jce_xfixed_trans, i + rjstart, &nc, &aj, NULL)); 1115 } 1116 1117 if (pdipm->nh) { 1118 /* Insert -grad h' */ 1119 PetscCall(MatGetRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, NULL)); 1120 PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality, &ranges)); 1121 proc = 0; 1122 for (j = 0; j < nc; j++) { 1123 /* find row ownership of */ 1124 while (aj[j] >= ranges[proc + 1]) proc++; 1125 nx_all = rranges[proc + 1] - rranges[proc]; 1126 col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 1127 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1128 } 1129 PetscCall(MatRestoreRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, NULL)); 1130 } 1131 1132 /* Insert Jci_xb^T' */ 1133 PetscCall(MatGetRow(Jci_xb_trans, i + rjstart, &nc, &aj, NULL)); 1134 PetscCall(MatGetOwnershipRanges(pdipm->Jci_xb, &ranges)); 1135 proc = 0; 1136 for (j = 0; j < nc; j++) { 1137 /* find row ownership of */ 1138 while (aj[j] >= ranges[proc + 1]) proc++; 1139 nx_all = rranges[proc + 1] - rranges[proc]; 1140 col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc]; 1141 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1142 } 1143 PetscCall(MatRestoreRow(Jci_xb_trans, i + rjstart, &nc, &aj, NULL)); 1144 } 1145 1146 /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */ 1147 if (pdipm->Ng) { 1148 PetscCall(MatGetOwnershipRange(tao->jacobian_equality, &rjstart, NULL)); 1149 for (i = 0; i < pdipm->ng; i++) { 1150 row = rstart + pdipm->off_lambdae + i; 1151 1152 PetscCall(MatGetRow(tao->jacobian_equality, i + rjstart, &nc, &aj, NULL)); 1153 proc = 0; 1154 for (j = 0; j < nc; j++) { 1155 while (aj[j] >= cranges[proc + 1]) proc++; 1156 col = aj[j] - cranges[proc] + Jranges[proc]; 1157 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); /* grad g */ 1158 } 1159 PetscCall(MatRestoreRow(tao->jacobian_equality, i + rjstart, &nc, &aj, NULL)); 1160 } 1161 } 1162 /* Jce_xfixed */ 1163 if (pdipm->Nxfixed) { 1164 PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, NULL)); 1165 for (i = 0; i < (pdipm->nce - pdipm->ng); i++) { 1166 row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1167 1168 PetscCall(MatGetRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, NULL)); 1169 PetscCheck(nc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "nc != 1"); 1170 1171 proc = 0; 1172 j = 0; 1173 while (cols[j] >= cranges[proc + 1]) proc++; 1174 col = cols[j] - cranges[proc] + Jranges[proc]; 1175 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1176 PetscCall(MatRestoreRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, NULL)); 1177 } 1178 } 1179 1180 /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */ 1181 if (pdipm->Nh) { 1182 PetscCall(MatGetOwnershipRange(tao->jacobian_inequality, &rjstart, NULL)); 1183 for (i = 0; i < pdipm->nh; i++) { 1184 row = rstart + pdipm->off_lambdai + i; 1185 1186 PetscCall(MatGetRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, NULL)); 1187 proc = 0; 1188 for (j = 0; j < nc; j++) { 1189 while (aj[j] >= cranges[proc + 1]) proc++; 1190 col = aj[j] - cranges[proc] + Jranges[proc]; 1191 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); /* grad h */ 1192 } 1193 PetscCall(MatRestoreRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, NULL)); 1194 } 1195 /* I */ 1196 for (i = 0; i < pdipm->nh; i++) { 1197 row = rstart + pdipm->off_lambdai + i; 1198 col = rstart + pdipm->off_z + i; 1199 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1200 } 1201 } 1202 1203 /* Jci_xb */ 1204 PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, NULL)); 1205 for (i = 0; i < (pdipm->nci - pdipm->nh); i++) { 1206 row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1207 1208 PetscCall(MatGetRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, NULL)); 1209 PetscCheck(nc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "nc != 1"); 1210 proc = 0; 1211 for (j = 0; j < nc; j++) { 1212 while (cols[j] >= cranges[proc + 1]) proc++; 1213 col = cols[j] - cranges[proc] + Jranges[proc]; 1214 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1215 } 1216 PetscCall(MatRestoreRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, NULL)); 1217 /* I */ 1218 col = rstart + pdipm->off_z + pdipm->nh + i; 1219 PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1220 } 1221 1222 /* 4-th Row block of KKT matrix: Z and Ci */ 1223 for (i = 0; i < pdipm->nci; i++) { 1224 row = rstart + pdipm->off_z + i; 1225 cols1[0] = rstart + pdipm->off_lambdai + i; 1226 cols1[1] = row; 1227 PetscCall(MatPreallocateSet(row, 2, cols1, dnz, onz)); 1228 } 1229 1230 /* diagonal entry */ 1231 for (i = 0; i < pdipm->n; i++) dnz[i]++; /* diagonal entry */ 1232 1233 /* Create KKT matrix */ 1234 PetscCall(MatCreate(comm, &J)); 1235 PetscCall(MatSetSizes(J, pdipm->n, pdipm->n, PETSC_DECIDE, PETSC_DECIDE)); 1236 PetscCall(MatSetFromOptions(J)); 1237 PetscCall(MatSeqAIJSetPreallocation(J, 0, dnz)); 1238 PetscCall(MatMPIAIJSetPreallocation(J, 0, dnz, 0, onz)); 1239 MatPreallocateEnd(dnz, onz); 1240 pdipm->K = J; 1241 1242 /* (8) Insert constant entries to K */ 1243 /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */ 1244 PetscCall(MatGetOwnershipRange(J, &rstart, &rend)); 1245 for (i = rstart; i < rend; i++) PetscCall(MatSetValue(J, i, i, 0.0, INSERT_VALUES)); 1246 /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */ 1247 if (pdipm->kkt_pd) { 1248 for (i = 0; i < pdipm->nh; i++) { 1249 row = rstart + i; 1250 PetscCall(MatSetValue(J, row, row, pdipm->deltaw, INSERT_VALUES)); 1251 } 1252 } 1253 1254 /* Row block of K: [ grad Ce, 0, 0, 0] */ 1255 if (pdipm->Nxfixed) { 1256 PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, NULL)); 1257 for (i = 0; i < (pdipm->nce - pdipm->ng); i++) { 1258 row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1259 1260 PetscCall(MatGetRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, &aa)); 1261 proc = 0; 1262 for (j = 0; j < nc; j++) { 1263 while (cols[j] >= cranges[proc + 1]) proc++; 1264 col = cols[j] - cranges[proc] + Jranges[proc]; 1265 PetscCall(MatSetValue(J, row, col, aa[j], INSERT_VALUES)); /* grad Ce */ 1266 PetscCall(MatSetValue(J, col, row, aa[j], INSERT_VALUES)); /* grad Ce' */ 1267 } 1268 PetscCall(MatRestoreRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, &aa)); 1269 } 1270 } 1271 1272 /* Row block of K: [ -grad Ci, 0, 0, I] */ 1273 PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, NULL)); 1274 for (i = 0; i < pdipm->nci - pdipm->nh; i++) { 1275 row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1276 1277 PetscCall(MatGetRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, &aa)); 1278 proc = 0; 1279 for (j = 0; j < nc; j++) { 1280 while (cols[j] >= cranges[proc + 1]) proc++; 1281 col = cols[j] - cranges[proc] + Jranges[proc]; 1282 PetscCall(MatSetValue(J, col, row, -aa[j], INSERT_VALUES)); 1283 PetscCall(MatSetValue(J, row, col, -aa[j], INSERT_VALUES)); 1284 } 1285 PetscCall(MatRestoreRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, &aa)); 1286 1287 col = rstart + pdipm->off_z + pdipm->nh + i; 1288 PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES)); 1289 } 1290 1291 for (i = 0; i < pdipm->nh; i++) { 1292 row = rstart + pdipm->off_lambdai + i; 1293 col = rstart + pdipm->off_z + i; 1294 PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES)); 1295 } 1296 1297 /* Row block of K: [ 0, 0, I, ...] */ 1298 for (i = 0; i < pdipm->nci; i++) { 1299 row = rstart + pdipm->off_z + i; 1300 col = rstart + pdipm->off_lambdai + i; 1301 PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES)); 1302 } 1303 1304 if (pdipm->Nxfixed) PetscCall(MatDestroy(&Jce_xfixed_trans)); 1305 PetscCall(MatDestroy(&Jci_xb_trans)); 1306 PetscCall(PetscFree3(ng_all, nh_all, Jranges)); 1307 1308 /* (9) Set up nonlinear solver SNES */ 1309 PetscCall(SNESSetFunction(pdipm->snes, NULL, TaoSNESFunction_PDIPM, (void *)tao)); 1310 PetscCall(SNESSetJacobian(pdipm->snes, J, J, TaoSNESJacobian_PDIPM, (void *)tao)); 1311 1312 if (pdipm->solve_reduced_kkt) { 1313 PC pc; 1314 PetscCall(KSPGetPC(tao->ksp, &pc)); 1315 PetscCall(PCSetType(pc, PCFIELDSPLIT)); 1316 PetscCall(PCFieldSplitSetType(pc, PC_COMPOSITE_SCHUR)); 1317 PetscCall(PCFieldSplitSetIS(pc, "2", pdipm->is2)); 1318 PetscCall(PCFieldSplitSetIS(pc, "1", pdipm->is1)); 1319 } 1320 PetscCall(SNESSetFromOptions(pdipm->snes)); 1321 1322 /* (10) Setup PCPostSetUp() for pdipm->solve_symmetric_kkt */ 1323 if (pdipm->solve_symmetric_kkt) { 1324 KSP ksp; 1325 PC pc; 1326 PetscBool isCHOL; 1327 1328 PetscCall(SNESGetKSP(pdipm->snes, &ksp)); 1329 PetscCall(KSPGetPC(ksp, &pc)); 1330 PetscCall(PCSetPostSetUp(pc, PCPostSetUp_PDIPM)); 1331 1332 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCCHOLESKY, &isCHOL)); 1333 if (isCHOL) { 1334 Mat Factor; 1335 1336 PetscCheck(PetscDefined(HAVE_MUMPS), PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Requires external package MUMPS"); 1337 PetscCall(PCFactorGetMatrix(pc, &Factor)); 1338 PetscCall(MatMumpsSetIcntl(Factor, 24, 1)); /* detection of null pivot rows */ 1339 if (size > 1) PetscCall(MatMumpsSetIcntl(Factor, 13, 1)); /* parallelism of the root node (enable ScaLAPACK) and its splitting */ 1340 } 1341 } 1342 PetscFunctionReturn(PETSC_SUCCESS); 1343 } 1344 1345 static PetscErrorCode TaoDestroy_PDIPM(Tao tao) 1346 { 1347 TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 1348 1349 PetscFunctionBegin; 1350 /* Freeing Vectors assocaiated with KKT (X) */ 1351 PetscCall(VecDestroy(&pdipm->x)); /* Solution x */ 1352 PetscCall(VecDestroy(&pdipm->lambdae)); /* Equality constraints lagrangian multiplier*/ 1353 PetscCall(VecDestroy(&pdipm->lambdai)); /* Inequality constraints lagrangian multiplier*/ 1354 PetscCall(VecDestroy(&pdipm->z)); /* Slack variables */ 1355 PetscCall(VecDestroy(&pdipm->X)); /* Big KKT system vector [x; lambdae; lambdai; z] */ 1356 1357 /* work vectors */ 1358 PetscCall(VecDestroy(&pdipm->lambdae_xfixed)); 1359 PetscCall(VecDestroy(&pdipm->lambdai_xb)); 1360 1361 /* Legrangian equality and inequality Vec */ 1362 PetscCall(VecDestroy(&pdipm->ce)); /* Vec of equality constraints */ 1363 PetscCall(VecDestroy(&pdipm->ci)); /* Vec of inequality constraints */ 1364 1365 /* Matrices */ 1366 PetscCall(MatDestroy(&pdipm->Jce_xfixed)); 1367 PetscCall(MatDestroy(&pdipm->Jci_xb)); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */ 1368 PetscCall(MatDestroy(&pdipm->K)); 1369 1370 /* Index Sets */ 1371 if (pdipm->Nxub) PetscCall(ISDestroy(&pdipm->isxub)); /* Finite upper bound only -inf < x < ub */ 1372 1373 if (pdipm->Nxlb) PetscCall(ISDestroy(&pdipm->isxlb)); /* Finite lower bound only lb <= x < inf */ 1374 1375 if (pdipm->Nxfixed) PetscCall(ISDestroy(&pdipm->isxfixed)); /* Fixed variables lb = x = ub */ 1376 1377 if (pdipm->Nxbox) PetscCall(ISDestroy(&pdipm->isxbox)); /* Boxed variables lb <= x <= ub */ 1378 1379 if (pdipm->Nxfree) PetscCall(ISDestroy(&pdipm->isxfree)); /* Free variables -inf <= x <= inf */ 1380 1381 if (pdipm->solve_reduced_kkt) { 1382 PetscCall(ISDestroy(&pdipm->is1)); 1383 PetscCall(ISDestroy(&pdipm->is2)); 1384 } 1385 1386 /* SNES */ 1387 PetscCall(SNESDestroy(&pdipm->snes)); /* Nonlinear solver */ 1388 PetscCall(PetscFree(pdipm->nce_all)); 1389 PetscCall(MatDestroy(&pdipm->jac_equality_trans)); 1390 PetscCall(MatDestroy(&pdipm->jac_inequality_trans)); 1391 1392 /* Destroy pdipm */ 1393 PetscCall(PetscFree(tao->data)); /* Holding locations of pdipm */ 1394 1395 /* Destroy Dual */ 1396 PetscCall(VecDestroy(&tao->DE)); /* equality dual */ 1397 PetscCall(VecDestroy(&tao->DI)); /* dinequality dual */ 1398 PetscFunctionReturn(PETSC_SUCCESS); 1399 } 1400 1401 static PetscErrorCode TaoSetFromOptions_PDIPM(Tao tao, PetscOptionItems PetscOptionsObject) 1402 { 1403 TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 1404 1405 PetscFunctionBegin; 1406 PetscOptionsHeadBegin(PetscOptionsObject, "PDIPM method for constrained optimization"); 1407 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)); 1408 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)); 1409 PetscCall(PetscOptionsBool("-tao_pdipm_solve_reduced_kkt", "Solve reduced KKT system using Schur-complement", NULL, pdipm->solve_reduced_kkt, &pdipm->solve_reduced_kkt, NULL)); 1410 PetscCall(PetscOptionsReal("-tao_pdipm_mu_update_factor", "Update scalar for barrier parameter (mu) update", NULL, pdipm->mu_update_factor, &pdipm->mu_update_factor, NULL)); 1411 PetscCall(PetscOptionsBool("-tao_pdipm_symmetric_kkt", "Solve non reduced symmetric KKT system", NULL, pdipm->solve_symmetric_kkt, &pdipm->solve_symmetric_kkt, NULL)); 1412 PetscCall(PetscOptionsBool("-tao_pdipm_kkt_shift_pd", "Add shifts to make KKT matrix positive definite", NULL, pdipm->kkt_pd, &pdipm->kkt_pd, NULL)); 1413 PetscOptionsHeadEnd(); 1414 PetscFunctionReturn(PETSC_SUCCESS); 1415 } 1416 1417 /*MC 1418 TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization. 1419 1420 Options Database Keys: 1421 + -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0) 1422 . -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0) 1423 . -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0) 1424 . -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system 1425 - -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite 1426 1427 Level: beginner 1428 1429 .seealso: `TAOPDIPM`, `Tao`, `TaoType` 1430 M*/ 1431 1432 PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao) 1433 { 1434 TAO_PDIPM *pdipm; 1435 PC pc; 1436 1437 PetscFunctionBegin; 1438 tao->ops->setup = TaoSetup_PDIPM; 1439 tao->ops->solve = TaoSolve_PDIPM; 1440 tao->ops->setfromoptions = TaoSetFromOptions_PDIPM; 1441 tao->ops->view = TaoView_PDIPM; 1442 tao->ops->destroy = TaoDestroy_PDIPM; 1443 1444 PetscCall(PetscNew(&pdipm)); 1445 tao->data = (void *)pdipm; 1446 1447 pdipm->nx = pdipm->Nx = 0; 1448 pdipm->nxfixed = pdipm->Nxfixed = 0; 1449 pdipm->nxlb = pdipm->Nxlb = 0; 1450 pdipm->nxub = pdipm->Nxub = 0; 1451 pdipm->nxbox = pdipm->Nxbox = 0; 1452 pdipm->nxfree = pdipm->Nxfree = 0; 1453 1454 pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0; 1455 pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0; 1456 pdipm->n = pdipm->N = 0; 1457 pdipm->mu = 1.0; 1458 pdipm->mu_update_factor = 0.1; 1459 1460 pdipm->deltaw = 0.0; 1461 pdipm->lastdeltaw = 3 * 1.e-4; 1462 pdipm->deltac = 0.0; 1463 pdipm->kkt_pd = PETSC_FALSE; 1464 1465 pdipm->push_init_slack = 1.0; 1466 pdipm->push_init_lambdai = 1.0; 1467 pdipm->solve_reduced_kkt = PETSC_FALSE; 1468 pdipm->solve_symmetric_kkt = PETSC_TRUE; 1469 1470 /* Override default settings (unless already changed) */ 1471 PetscCall(TaoParametersInitialize(tao)); 1472 PetscObjectParameterSetDefault(tao, max_it, 200); 1473 PetscObjectParameterSetDefault(tao, max_funcs, 500); 1474 1475 PetscCall(SNESCreate(((PetscObject)tao)->comm, &pdipm->snes)); 1476 PetscCall(SNESSetOptionsPrefix(pdipm->snes, tao->hdr.prefix)); 1477 PetscCall(SNESGetKSP(pdipm->snes, &tao->ksp)); 1478 PetscCall(PetscObjectReference((PetscObject)tao->ksp)); 1479 PetscCall(KSPGetPC(tao->ksp, &pc)); 1480 PetscCall(PCSetApplicationContext(pc, (void *)tao)); 1481 PetscFunctionReturn(PETSC_SUCCESS); 1482 } 1483