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