1ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/ 2c4b75bccSAlp Dener #include <petsc/private/vecimpl.h> 3af0996ceSBarry Smith #include <petsc/private/taoimpl.h> 4aaa7dc30SBarry Smith #include <../src/tao/matrix/submatfree.h> 5a7e14dcfSSatish Balay 6b98f30f2SJason Sarich /*@C 7b98f30f2SJason Sarich TaoVecGetSubVec - Gets a subvector using the IS 8a7e14dcfSSatish Balay 9a7e14dcfSSatish Balay Input Parameters: 10a7e14dcfSSatish Balay + vfull - the full matrix 11a7e14dcfSSatish Balay . is - the index set for the subvector 12a7e14dcfSSatish Balay . reduced_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK, TAO_SUBSET_MATRIXFREE) 13a7e14dcfSSatish Balay - maskvalue - the value to set the unused vector elements to (for TAO_SUBSET_MASK or TAO_SUBSET_MATRIXFREE) 14a7e14dcfSSatish Balay 1597bb3fdcSJose E. Roman Output Parameter: 16a7e14dcfSSatish Balay . vreduced - the subvector 17a7e14dcfSSatish Balay 181eb8069cSJason Sarich Notes: 19a7e14dcfSSatish Balay maskvalue should usually be 0.0, unless a pointwise divide will be used. 201eb8069cSJason Sarich 21b6a6cedcSAlp Dener Level: developer 22a7e14dcfSSatish Balay @*/ 239371c9d4SSatish Balay PetscErrorCode TaoVecGetSubVec(Vec vfull, IS is, TaoSubsetType reduced_type, PetscReal maskvalue, Vec *vreduced) { 24a7e14dcfSSatish Balay PetscInt nfull, nreduced, nreduced_local, rlow, rhigh, flow, fhigh; 25a7e14dcfSSatish Balay PetscInt i, nlocal; 26a7e14dcfSSatish Balay PetscReal *fv, *rv; 27a7e14dcfSSatish Balay const PetscInt *s; 28a7e14dcfSSatish Balay IS ident; 29a7e14dcfSSatish Balay VecType vtype; 30a7e14dcfSSatish Balay VecScatter scatter; 31a7e14dcfSSatish Balay MPI_Comm comm; 32a7e14dcfSSatish Balay 33a7e14dcfSSatish Balay PetscFunctionBegin; 34a7e14dcfSSatish Balay PetscValidHeaderSpecific(vfull, VEC_CLASSID, 1); 35a7e14dcfSSatish Balay PetscValidHeaderSpecific(is, IS_CLASSID, 2); 36a7e14dcfSSatish Balay 379566063dSJacob Faibussowitsch PetscCall(VecGetSize(vfull, &nfull)); 389566063dSJacob Faibussowitsch PetscCall(ISGetSize(is, &nreduced)); 39a7e14dcfSSatish Balay 40a7e14dcfSSatish Balay if (nreduced == nfull) { 419566063dSJacob Faibussowitsch PetscCall(VecDestroy(vreduced)); 429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vfull, vreduced)); 439566063dSJacob Faibussowitsch PetscCall(VecCopy(vfull, *vreduced)); 44a7e14dcfSSatish Balay } else { 45a7e14dcfSSatish Balay switch (reduced_type) { 46a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 479566063dSJacob Faibussowitsch PetscCall(VecGetType(vfull, &vtype)); 489566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vfull, &flow, &fhigh)); 499566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &nreduced_local)); 509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)vfull, &comm)); 511baa6e33SBarry Smith if (*vreduced) PetscCall(VecDestroy(vreduced)); 529566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, vreduced)); 539566063dSJacob Faibussowitsch PetscCall(VecSetType(*vreduced, vtype)); 54a7e14dcfSSatish Balay 559566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*vreduced, nreduced_local, nreduced)); 569566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(*vreduced, &rlow, &rhigh)); 579566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, nreduced_local, rlow, 1, &ident)); 589566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vfull, is, *vreduced, ident, &scatter)); 599566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter, vfull, *vreduced, INSERT_VALUES, SCATTER_FORWARD)); 609566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter, vfull, *vreduced, INSERT_VALUES, SCATTER_FORWARD)); 619566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scatter)); 629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ident)); 63a7e14dcfSSatish Balay break; 64a7e14dcfSSatish Balay 65a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 66a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 67a7e14dcfSSatish Balay /* vr[i] = vf[i] if i in is 68a7e14dcfSSatish Balay vr[i] = 0 otherwise */ 69*48a46eb9SPierre Jolivet if (!*vreduced) PetscCall(VecDuplicate(vfull, vreduced)); 70a7e14dcfSSatish Balay 719566063dSJacob Faibussowitsch PetscCall(VecSet(*vreduced, maskvalue)); 729566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &nlocal)); 739566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vfull, &flow, &fhigh)); 749566063dSJacob Faibussowitsch PetscCall(VecGetArray(vfull, &fv)); 759566063dSJacob Faibussowitsch PetscCall(VecGetArray(*vreduced, &rv)); 769566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &s)); 7763a3b9bcSJacob Faibussowitsch PetscCheck(nlocal <= (fhigh - flow), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "IS local size %" PetscInt_FMT " > Vec local size %" PetscInt_FMT, nlocal, fhigh - flow); 789371c9d4SSatish Balay for (i = 0; i < nlocal; ++i) { rv[s[i] - flow] = fv[s[i] - flow]; } 799566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &s)); 809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vfull, &fv)); 819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*vreduced, &rv)); 82a7e14dcfSSatish Balay break; 83a7e14dcfSSatish Balay } 84a7e14dcfSSatish Balay } 85a7e14dcfSSatish Balay PetscFunctionReturn(0); 86a7e14dcfSSatish Balay } 87a7e14dcfSSatish Balay 88b98f30f2SJason Sarich /*@C 89b98f30f2SJason Sarich TaoMatGetSubMat - Gets a submatrix using the IS 90a7e14dcfSSatish Balay 91a7e14dcfSSatish Balay Input Parameters: 92a7e14dcfSSatish Balay + M - the full matrix (n x n) 93a7e14dcfSSatish Balay . is - the index set for the submatrix (both row and column index sets need to be the same) 94a7e14dcfSSatish Balay . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option 95b6a6cedcSAlp Dener - subset_type <TAO_SUBSET_SUBVEC,TAO_SUBSET_MASK,TAO_SUBSET_MATRIXFREE> - the method TAO is using for subsetting 96a7e14dcfSSatish Balay 9797bb3fdcSJose E. Roman Output Parameter: 98a7e14dcfSSatish Balay . Msub - the submatrix 99b6a6cedcSAlp Dener 100b6a6cedcSAlp Dener Level: developer 101a7e14dcfSSatish Balay @*/ 1029371c9d4SSatish Balay PetscErrorCode TaoMatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub) { 103a7e14dcfSSatish Balay IS iscomp; 10462675beeSAlp Dener PetscBool flg = PETSC_TRUE; 10553506e15SBarry Smith 106a7e14dcfSSatish Balay PetscFunctionBegin; 107a7e14dcfSSatish Balay PetscValidHeaderSpecific(M, MAT_CLASSID, 1); 108a7e14dcfSSatish Balay PetscValidHeaderSpecific(is, IS_CLASSID, 2); 1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(Msub)); 110a7e14dcfSSatish Balay switch (subset_type) { 1119371c9d4SSatish Balay case TAO_SUBSET_SUBVEC: PetscCall(MatCreateSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub)); break; 112a7e14dcfSSatish Balay 113a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 114a7e14dcfSSatish Balay /* Get Reduced Hessian 115a7e14dcfSSatish Balay Msub[i,j] = M[i,j] if i,j in Free_Local or i==j 116a7e14dcfSSatish Balay Msub[i,j] = 0 if i!=j and i or j not in Free_Local 117a7e14dcfSSatish Balay */ 118d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)M); 1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-overwrite_hessian", "modify the existing hessian matrix when computing submatrices", "TaoSubsetType", flg, &flg, NULL)); 120d0609cedSBarry Smith PetscOptionsEnd(); 1218afaa268SBarry Smith if (flg) { 1229566063dSJacob Faibussowitsch PetscCall(MatDuplicate(M, MAT_COPY_VALUES, Msub)); 123a7e14dcfSSatish Balay } else { 124a7e14dcfSSatish Balay /* Act on hessian directly (default) */ 1259566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)M)); 126a7e14dcfSSatish Balay *Msub = M; 127a7e14dcfSSatish Balay } 128a7e14dcfSSatish Balay /* Save the diagonal to temporary vector */ 1299566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(*Msub, v1)); 130a7e14dcfSSatish Balay 131a7e14dcfSSatish Balay /* Zero out rows and columns */ 1329566063dSJacob Faibussowitsch PetscCall(ISComplementVec(is, v1, &iscomp)); 133a7e14dcfSSatish Balay 134a7e14dcfSSatish Balay /* Use v1 instead of 0 here because of PETSc bug */ 1359566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsIS(*Msub, iscomp, 1.0, v1, v1)); 136a7e14dcfSSatish Balay 1379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscomp)); 138a7e14dcfSSatish Balay break; 139a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 1409566063dSJacob Faibussowitsch PetscCall(ISComplementVec(is, v1, &iscomp)); 1419566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrixFree(M, iscomp, iscomp, Msub)); 1429566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscomp)); 143a7e14dcfSSatish Balay break; 144a7e14dcfSSatish Balay } 145a7e14dcfSSatish Balay PetscFunctionReturn(0); 146a7e14dcfSSatish Balay } 1472f75a4aaSAlp Dener 1482f75a4aaSAlp Dener /*@C 1492f75a4aaSAlp Dener TaoEstimateActiveBounds - Generates index sets for variables at the lower and upper 1502f75a4aaSAlp Dener bounds, as well as fixed variables where lower and upper bounds equal each other. 1512f75a4aaSAlp Dener 1522f75a4aaSAlp Dener Input Parameters: 1532f75a4aaSAlp Dener + X - solution vector 1542f75a4aaSAlp Dener . XL - lower bound vector 1552f75a4aaSAlp Dener . XU - upper bound vector 1562f75a4aaSAlp Dener . G - unprojected gradient 1570a4511e9SAlp Dener . S - step direction with which the active bounds will be estimated 158b6a6cedcSAlp Dener . W - work vector of type and size of X 1590a4511e9SAlp Dener - steplen - the step length at which the active bounds will be estimated (needs to be conservative) 1602f75a4aaSAlp Dener 1612f75a4aaSAlp Dener Output Parameters: 16276be6f4fSStefano Zampini + bound_tol - tolerance for the bound estimation 1632f75a4aaSAlp Dener . active_lower - index set for active variables at the lower bound 1642f75a4aaSAlp Dener . active_upper - index set for active variables at the upper bound 1652f75a4aaSAlp Dener . active_fixed - index set for fixed variables 1662f75a4aaSAlp Dener . active - index set for all active variables 167b6a6cedcSAlp Dener - inactive - complementary index set for inactive variables 1680a4511e9SAlp Dener 1690a4511e9SAlp Dener Notes: 1700a4511e9SAlp Dener This estimation is based on Bertsekas' method, with a built in diagonal scaling value of 1.0e-3. 1710a4511e9SAlp Dener 172b6a6cedcSAlp Dener Level: developer 1732f75a4aaSAlp Dener @*/ 1749371c9d4SSatish Balay PetscErrorCode TaoEstimateActiveBounds(Vec X, Vec XL, Vec XU, Vec G, Vec S, Vec W, PetscReal steplen, PetscReal *bound_tol, IS *active_lower, IS *active_upper, IS *active_fixed, IS *active, IS *inactive) { 1752f75a4aaSAlp Dener PetscReal wnorm; 17689da521bSAlp Dener PetscReal zero = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0); 177c4b75bccSAlp Dener PetscInt i, n_isl = 0, n_isu = 0, n_isf = 0, n_isa = 0, n_isi = 0; 178c4b75bccSAlp Dener PetscInt N_isl, N_isu, N_isf, N_isa, N_isi; 179c4b75bccSAlp Dener PetscInt n, low, high, nDiff; 180c4b75bccSAlp Dener PetscInt *isl = NULL, *isu = NULL, *isf = NULL, *isa = NULL, *isi = NULL; 1812f75a4aaSAlp Dener const PetscScalar *xl, *xu, *x, *g; 1826519dc0cSAlp Dener MPI_Comm comm = PetscObjectComm((PetscObject)X); 1832f75a4aaSAlp Dener 1842f75a4aaSAlp Dener PetscFunctionBegin; 185c4b75bccSAlp Dener PetscValidHeaderSpecific(X, VEC_CLASSID, 1); 18676be6f4fSStefano Zampini if (XL) PetscValidHeaderSpecific(XL, VEC_CLASSID, 2); 18776be6f4fSStefano Zampini if (XU) PetscValidHeaderSpecific(XU, VEC_CLASSID, 3); 188c4b75bccSAlp Dener PetscValidHeaderSpecific(G, VEC_CLASSID, 4); 189c4b75bccSAlp Dener PetscValidHeaderSpecific(S, VEC_CLASSID, 5); 190c4b75bccSAlp Dener PetscValidHeaderSpecific(W, VEC_CLASSID, 6); 191c4b75bccSAlp Dener 19276be6f4fSStefano Zampini if (XL) PetscCheckSameType(X, 1, XL, 2); 19376be6f4fSStefano Zampini if (XU) PetscCheckSameType(X, 1, XU, 3); 194c4b75bccSAlp Dener PetscCheckSameType(X, 1, G, 4); 195c4b75bccSAlp Dener PetscCheckSameType(X, 1, S, 5); 196c4b75bccSAlp Dener PetscCheckSameType(X, 1, W, 6); 19776be6f4fSStefano Zampini if (XL) PetscCheckSameComm(X, 1, XL, 2); 19876be6f4fSStefano Zampini if (XU) PetscCheckSameComm(X, 1, XU, 3); 199c4b75bccSAlp Dener PetscCheckSameComm(X, 1, G, 4); 200c4b75bccSAlp Dener PetscCheckSameComm(X, 1, S, 5); 201c4b75bccSAlp Dener PetscCheckSameComm(X, 1, W, 6); 20276be6f4fSStefano Zampini if (XL) VecCheckSameSize(X, 1, XL, 2); 20376be6f4fSStefano Zampini if (XU) VecCheckSameSize(X, 1, XU, 3); 204c4b75bccSAlp Dener VecCheckSameSize(X, 1, G, 4); 205c4b75bccSAlp Dener VecCheckSameSize(X, 1, S, 5); 206c4b75bccSAlp Dener VecCheckSameSize(X, 1, W, 6); 207c4b75bccSAlp Dener 2082f75a4aaSAlp Dener /* Update the tolerance for bound detection (this is based on Bertsekas' method) */ 2099566063dSJacob Faibussowitsch PetscCall(VecCopy(X, W)); 2109566063dSJacob Faibussowitsch PetscCall(VecAXPBY(W, steplen, 1.0, S)); 2119566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(W, XL, XU, 0.0, &nDiff, W)); 2129566063dSJacob Faibussowitsch PetscCall(VecAXPBY(W, 1.0, -1.0, X)); 2139566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 2142f75a4aaSAlp Dener *bound_tol = PetscMin(*bound_tol, wnorm); 2152f75a4aaSAlp Dener 21676be6f4fSStefano Zampini /* Clear all index sets */ 21776be6f4fSStefano Zampini PetscCall(ISDestroy(active_lower)); 21876be6f4fSStefano Zampini PetscCall(ISDestroy(active_upper)); 21976be6f4fSStefano Zampini PetscCall(ISDestroy(active_fixed)); 22076be6f4fSStefano Zampini PetscCall(ISDestroy(active)); 22176be6f4fSStefano Zampini PetscCall(ISDestroy(inactive)); 22276be6f4fSStefano Zampini 2239566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &low, &high)); 2249566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 22576be6f4fSStefano Zampini if (!XL && !XU) { 22676be6f4fSStefano Zampini PetscCall(ISCreateStride(comm, n, low, 1, inactive)); 22776be6f4fSStefano Zampini PetscFunctionReturn(0); 22876be6f4fSStefano Zampini } 2292f75a4aaSAlp Dener if (n > 0) { 2309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XL, &xl)); 2329566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XU, &xu)); 2339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(G, &g)); 2342f75a4aaSAlp Dener 2352f75a4aaSAlp Dener /* Loop over variables and categorize the indexes */ 2369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isl)); 2379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isu)); 2389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isf)); 2399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isa)); 2409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isi)); 241c4b75bccSAlp Dener for (i = 0; i < n; ++i) { 24266b4d56fSTodd Munson if (xl[i] == xu[i]) { 243c4b75bccSAlp Dener /* Fixed variables */ 2449371c9d4SSatish Balay isf[n_isf] = low + i; 2459371c9d4SSatish Balay ++n_isf; 2469371c9d4SSatish Balay isa[n_isa] = low + i; 2479371c9d4SSatish Balay ++n_isa; 24876be6f4fSStefano Zampini } else if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + *bound_tol && g[i] > zero) { 249c4b75bccSAlp Dener /* Lower bounded variables */ 2509371c9d4SSatish Balay isl[n_isl] = low + i; 2519371c9d4SSatish Balay ++n_isl; 2529371c9d4SSatish Balay isa[n_isa] = low + i; 2539371c9d4SSatish Balay ++n_isa; 25476be6f4fSStefano Zampini } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - *bound_tol && g[i] < zero) { 255c4b75bccSAlp Dener /* Upper bounded variables */ 2569371c9d4SSatish Balay isu[n_isu] = low + i; 2579371c9d4SSatish Balay ++n_isu; 2589371c9d4SSatish Balay isa[n_isa] = low + i; 2599371c9d4SSatish Balay ++n_isa; 260c4b75bccSAlp Dener } else { 261c4b75bccSAlp Dener /* Inactive variables */ 2629371c9d4SSatish Balay isi[n_isi] = low + i; 2639371c9d4SSatish Balay ++n_isi; 2642f75a4aaSAlp Dener } 2652f75a4aaSAlp Dener } 2662f75a4aaSAlp Dener 2679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XL, &xl)); 2699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XU, &xu)); 2709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(G, &g)); 2712f75a4aaSAlp Dener } 2722f75a4aaSAlp Dener 273c4b75bccSAlp Dener /* Collect global sizes */ 2741c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&n_isl, &N_isl, 1, MPIU_INT, MPI_SUM, comm)); 2751c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&n_isu, &N_isu, 1, MPIU_INT, MPI_SUM, comm)); 2761c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&n_isf, &N_isf, 1, MPIU_INT, MPI_SUM, comm)); 2771c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&n_isa, &N_isa, 1, MPIU_INT, MPI_SUM, comm)); 2781c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&n_isi, &N_isi, 1, MPIU_INT, MPI_SUM, comm)); 279c4b75bccSAlp Dener 280c4b75bccSAlp Dener /* Create index set for lower bounded variables */ 281c4b75bccSAlp Dener if (N_isl > 0) { 2829566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_isl, isl, PETSC_OWN_POINTER, active_lower)); 2837e9273c4SAlp Dener } else { 2849566063dSJacob Faibussowitsch PetscCall(PetscFree(isl)); 285c4b75bccSAlp Dener } 286c4b75bccSAlp Dener /* Create index set for upper bounded variables */ 287c4b75bccSAlp Dener if (N_isu > 0) { 2889566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_isu, isu, PETSC_OWN_POINTER, active_upper)); 2897e9273c4SAlp Dener } else { 2909566063dSJacob Faibussowitsch PetscCall(PetscFree(isu)); 291c4b75bccSAlp Dener } 292c4b75bccSAlp Dener /* Create index set for fixed variables */ 293c4b75bccSAlp Dener if (N_isf > 0) { 2949566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_isf, isf, PETSC_OWN_POINTER, active_fixed)); 2957e9273c4SAlp Dener } else { 2969566063dSJacob Faibussowitsch PetscCall(PetscFree(isf)); 297c4b75bccSAlp Dener } 298c4b75bccSAlp Dener /* Create index set for all actively bounded variables */ 299c4b75bccSAlp Dener if (N_isa > 0) { 3009566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_isa, isa, PETSC_OWN_POINTER, active)); 3017e9273c4SAlp Dener } else { 3029566063dSJacob Faibussowitsch PetscCall(PetscFree(isa)); 303c4b75bccSAlp Dener } 304c4b75bccSAlp Dener /* Create index set for all inactive variables */ 305c4b75bccSAlp Dener if (N_isi > 0) { 3069566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_isi, isi, PETSC_OWN_POINTER, inactive)); 3077e9273c4SAlp Dener } else { 3089566063dSJacob Faibussowitsch PetscCall(PetscFree(isi)); 309c4b75bccSAlp Dener } 3102f75a4aaSAlp Dener PetscFunctionReturn(0); 3112f75a4aaSAlp Dener } 3122f75a4aaSAlp Dener 3132f75a4aaSAlp Dener /*@C 3142f75a4aaSAlp Dener TaoBoundStep - Ensures the correct zero or adjusted step direction 3152f75a4aaSAlp Dener values for active variables. 3162f75a4aaSAlp Dener 3172f75a4aaSAlp Dener Input Parameters: 3182f75a4aaSAlp Dener + X - solution vector 3192f75a4aaSAlp Dener . XL - lower bound vector 3202f75a4aaSAlp Dener . XU - upper bound vector 3212f75a4aaSAlp Dener . active_lower - index set for lower bounded active variables 3222f75a4aaSAlp Dener . active_upper - index set for lower bounded active variables 323b6a6cedcSAlp Dener . active_fixed - index set for fixed active variables 324b6a6cedcSAlp Dener - scale - amplification factor for the step that needs to be taken on actively bounded variables 3252f75a4aaSAlp Dener 32697bb3fdcSJose E. Roman Output Parameter: 3272f75a4aaSAlp Dener . S - step direction to be modified 328b6a6cedcSAlp Dener 329b6a6cedcSAlp Dener Level: developer 3302f75a4aaSAlp Dener @*/ 3319371c9d4SSatish Balay PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, PetscReal scale, Vec S) { 3322f75a4aaSAlp Dener Vec step_lower, step_upper, step_fixed; 3332f75a4aaSAlp Dener Vec x_lower, x_upper; 3342f75a4aaSAlp Dener Vec bound_lower, bound_upper; 3352f75a4aaSAlp Dener 3362f75a4aaSAlp Dener PetscFunctionBegin; 3372f75a4aaSAlp Dener /* Adjust step for variables at the estimated lower bound */ 3382f75a4aaSAlp Dener if (active_lower) { 3399566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(S, active_lower, &step_lower)); 3409566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, active_lower, &x_lower)); 3419566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(XL, active_lower, &bound_lower)); 3429566063dSJacob Faibussowitsch PetscCall(VecCopy(bound_lower, step_lower)); 3439566063dSJacob Faibussowitsch PetscCall(VecAXPY(step_lower, -1.0, x_lower)); 3449566063dSJacob Faibussowitsch PetscCall(VecScale(step_lower, scale)); 3459566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(S, active_lower, &step_lower)); 3469566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, active_lower, &x_lower)); 3479566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(XL, active_lower, &bound_lower)); 3482f75a4aaSAlp Dener } 3492f75a4aaSAlp Dener 3502f75a4aaSAlp Dener /* Adjust step for the variables at the estimated upper bound */ 3512f75a4aaSAlp Dener if (active_upper) { 3529566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(S, active_upper, &step_upper)); 3539566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, active_upper, &x_upper)); 3549566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(XU, active_upper, &bound_upper)); 3559566063dSJacob Faibussowitsch PetscCall(VecCopy(bound_upper, step_upper)); 3569566063dSJacob Faibussowitsch PetscCall(VecAXPY(step_upper, -1.0, x_upper)); 3579566063dSJacob Faibussowitsch PetscCall(VecScale(step_upper, scale)); 3589566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(S, active_upper, &step_upper)); 3599566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, active_upper, &x_upper)); 3609566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(XU, active_upper, &bound_upper)); 3612f75a4aaSAlp Dener } 3622f75a4aaSAlp Dener 3632f75a4aaSAlp Dener /* Zero out step for fixed variables */ 3642f75a4aaSAlp Dener if (active_fixed) { 3659566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(S, active_fixed, &step_fixed)); 3669566063dSJacob Faibussowitsch PetscCall(VecSet(step_fixed, 0.0)); 3679566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(S, active_fixed, &step_fixed)); 3682f75a4aaSAlp Dener } 3692f75a4aaSAlp Dener PetscFunctionReturn(0); 3702f75a4aaSAlp Dener } 371c4b75bccSAlp Dener 372c4b75bccSAlp Dener /*@C 373b6a6cedcSAlp Dener TaoBoundSolution - Ensures that the solution vector is snapped into the bounds within a given tolerance. 374c4b75bccSAlp Dener 375f0fc11ceSJed Brown Collective on Vec 376f0fc11ceSJed Brown 377c4b75bccSAlp Dener Input Parameters: 378b6a6cedcSAlp Dener + X - solution vector 379b6a6cedcSAlp Dener . XL - lower bound vector 380c4b75bccSAlp Dener . XU - upper bound vector 381f0fc11ceSJed Brown - bound_tol - absolute tolerance in enforcing the bound 382c4b75bccSAlp Dener 383c4b75bccSAlp Dener Output Parameters: 384f0fc11ceSJed Brown + nDiff - total number of vector entries that have been bounded 385f0fc11ceSJed Brown - Xout - modified solution vector satisfying bounds to bound_tol 386b6a6cedcSAlp Dener 387b6a6cedcSAlp Dener Level: developer 388f0fc11ceSJed Brown 389db781477SPatrick Sanan .seealso: `TAOBNCG`, `TAOBNTL`, `TAOBNTR` 390c4b75bccSAlp Dener @*/ 3919371c9d4SSatish Balay PetscErrorCode TaoBoundSolution(Vec X, Vec XL, Vec XU, PetscReal bound_tol, PetscInt *nDiff, Vec Xout) { 392c4b75bccSAlp Dener PetscInt i, n, low, high, nDiff_loc = 0; 3933b063059SAlp Dener PetscScalar *xout; 3943b063059SAlp Dener const PetscScalar *x, *xl, *xu; 395c4b75bccSAlp Dener 396c4b75bccSAlp Dener PetscFunctionBegin; 397c4b75bccSAlp Dener PetscValidHeaderSpecific(X, VEC_CLASSID, 1); 39876be6f4fSStefano Zampini if (XL) PetscValidHeaderSpecific(XL, VEC_CLASSID, 2); 39976be6f4fSStefano Zampini if (XU) PetscValidHeaderSpecific(XU, VEC_CLASSID, 3); 400064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(Xout, VEC_CLASSID, 6); 40176be6f4fSStefano Zampini if (!XL && !XU) { 40276be6f4fSStefano Zampini PetscCall(VecCopy(X, Xout)); 40376be6f4fSStefano Zampini *nDiff = 0.0; 40476be6f4fSStefano Zampini PetscFunctionReturn(0); 40576be6f4fSStefano Zampini } 406c4b75bccSAlp Dener PetscCheckSameType(X, 1, XL, 2); 407c4b75bccSAlp Dener PetscCheckSameType(X, 1, XU, 3); 408064a246eSJacob Faibussowitsch PetscCheckSameType(X, 1, Xout, 6); 409c4b75bccSAlp Dener PetscCheckSameComm(X, 1, XL, 2); 410c4b75bccSAlp Dener PetscCheckSameComm(X, 1, XU, 3); 411064a246eSJacob Faibussowitsch PetscCheckSameComm(X, 1, Xout, 6); 412c4b75bccSAlp Dener VecCheckSameSize(X, 1, XL, 2); 413c4b75bccSAlp Dener VecCheckSameSize(X, 1, XU, 3); 4143b063059SAlp Dener VecCheckSameSize(X, 1, Xout, 4); 415c4b75bccSAlp Dener 4169566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &low, &high)); 4179566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 418c4b75bccSAlp Dener if (n > 0) { 4199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 4209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XL, &xl)); 4219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XU, &xu)); 4229566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xout, &xout)); 423c4b75bccSAlp Dener 424c4b75bccSAlp Dener for (i = 0; i < n; ++i) { 42576be6f4fSStefano Zampini if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + bound_tol) { 4269371c9d4SSatish Balay xout[i] = xl[i]; 4279371c9d4SSatish Balay ++nDiff_loc; 42876be6f4fSStefano Zampini } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - bound_tol) { 4299371c9d4SSatish Balay xout[i] = xu[i]; 4309371c9d4SSatish Balay ++nDiff_loc; 431c4b75bccSAlp Dener } 432c4b75bccSAlp Dener } 433c4b75bccSAlp Dener 4349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 4359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XL, &xl)); 4369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XU, &xu)); 4379566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xout, &xout)); 438c4b75bccSAlp Dener } 4391c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&nDiff_loc, nDiff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)X))); 440c4b75bccSAlp Dener PetscFunctionReturn(0); 441c4b75bccSAlp Dener } 442