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