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