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