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 @*/ 23d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoVecGetSubVec(Vec vfull, IS is, TaoSubsetType reduced_type, PetscReal maskvalue, Vec *vreduced) 24d71ae5a4SJacob Faibussowitsch { 25a7e14dcfSSatish Balay PetscInt nfull, nreduced, nreduced_local, rlow, rhigh, flow, fhigh; 26a7e14dcfSSatish Balay PetscInt i, nlocal; 27a7e14dcfSSatish Balay PetscReal *fv, *rv; 28a7e14dcfSSatish Balay const PetscInt *s; 29a7e14dcfSSatish Balay IS ident; 30a7e14dcfSSatish Balay VecType vtype; 31a7e14dcfSSatish Balay VecScatter scatter; 32a7e14dcfSSatish Balay MPI_Comm comm; 33a7e14dcfSSatish Balay 34a7e14dcfSSatish Balay PetscFunctionBegin; 35a7e14dcfSSatish Balay PetscValidHeaderSpecific(vfull, VEC_CLASSID, 1); 36a7e14dcfSSatish Balay PetscValidHeaderSpecific(is, IS_CLASSID, 2); 37a7e14dcfSSatish Balay 389566063dSJacob Faibussowitsch PetscCall(VecGetSize(vfull, &nfull)); 399566063dSJacob Faibussowitsch PetscCall(ISGetSize(is, &nreduced)); 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay if (nreduced == nfull) { 429566063dSJacob Faibussowitsch PetscCall(VecDestroy(vreduced)); 439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vfull, vreduced)); 449566063dSJacob Faibussowitsch PetscCall(VecCopy(vfull, *vreduced)); 45a7e14dcfSSatish Balay } else { 46a7e14dcfSSatish Balay switch (reduced_type) { 47a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 489566063dSJacob Faibussowitsch PetscCall(VecGetType(vfull, &vtype)); 499566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vfull, &flow, &fhigh)); 509566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &nreduced_local)); 519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)vfull, &comm)); 521baa6e33SBarry Smith if (*vreduced) PetscCall(VecDestroy(vreduced)); 539566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, vreduced)); 549566063dSJacob Faibussowitsch PetscCall(VecSetType(*vreduced, vtype)); 55a7e14dcfSSatish Balay 569566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*vreduced, nreduced_local, nreduced)); 579566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(*vreduced, &rlow, &rhigh)); 589566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, nreduced_local, rlow, 1, &ident)); 599566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vfull, is, *vreduced, ident, &scatter)); 609566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter, vfull, *vreduced, INSERT_VALUES, SCATTER_FORWARD)); 619566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter, vfull, *vreduced, INSERT_VALUES, SCATTER_FORWARD)); 629566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scatter)); 639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ident)); 64a7e14dcfSSatish Balay break; 65a7e14dcfSSatish Balay 66a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 67a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 68a7e14dcfSSatish Balay /* vr[i] = vf[i] if i in is 69a7e14dcfSSatish Balay vr[i] = 0 otherwise */ 7048a46eb9SPierre Jolivet if (!*vreduced) PetscCall(VecDuplicate(vfull, vreduced)); 71a7e14dcfSSatish Balay 729566063dSJacob Faibussowitsch PetscCall(VecSet(*vreduced, maskvalue)); 739566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &nlocal)); 749566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vfull, &flow, &fhigh)); 759566063dSJacob Faibussowitsch PetscCall(VecGetArray(vfull, &fv)); 769566063dSJacob Faibussowitsch PetscCall(VecGetArray(*vreduced, &rv)); 779566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &s)); 7863a3b9bcSJacob 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); 79ad540459SPierre Jolivet for (i = 0; i < nlocal; ++i) rv[s[i] - flow] = fv[s[i] - flow]; 809566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &s)); 819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vfull, &fv)); 829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*vreduced, &rv)); 83a7e14dcfSSatish Balay break; 84a7e14dcfSSatish Balay } 85a7e14dcfSSatish Balay } 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87a7e14dcfSSatish Balay } 88a7e14dcfSSatish Balay 89b98f30f2SJason Sarich /*@C 90*2fe279fdSBarry Smith TaoMatGetSubMat - Gets a submatrix using the `IS` 91a7e14dcfSSatish Balay 92a7e14dcfSSatish Balay Input Parameters: 93a7e14dcfSSatish Balay + M - the full matrix (n x n) 94a7e14dcfSSatish Balay . is - the index set for the submatrix (both row and column index sets need to be the same) 95a7e14dcfSSatish Balay . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option 96*2fe279fdSBarry Smith - subset_type <`TAO_SUBSET_SUBVEC`, `TAO_SUBSET_MASK`, `TAO_SUBSET_MATRIXFREE`> - the method `Tao` is using for subsetting 97a7e14dcfSSatish Balay 9897bb3fdcSJose E. Roman Output Parameter: 99a7e14dcfSSatish Balay . Msub - the submatrix 100b6a6cedcSAlp Dener 101b6a6cedcSAlp Dener Level: developer 102a7e14dcfSSatish Balay @*/ 103d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoMatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub) 104d71ae5a4SJacob Faibussowitsch { 105a7e14dcfSSatish Balay IS iscomp; 10662675beeSAlp Dener PetscBool flg = PETSC_TRUE; 10753506e15SBarry Smith 108a7e14dcfSSatish Balay PetscFunctionBegin; 109a7e14dcfSSatish Balay PetscValidHeaderSpecific(M, MAT_CLASSID, 1); 110a7e14dcfSSatish Balay PetscValidHeaderSpecific(is, IS_CLASSID, 2); 1119566063dSJacob Faibussowitsch PetscCall(MatDestroy(Msub)); 112a7e14dcfSSatish Balay switch (subset_type) { 113d71ae5a4SJacob Faibussowitsch case TAO_SUBSET_SUBVEC: 114d71ae5a4SJacob Faibussowitsch PetscCall(MatCreateSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub)); 115d71ae5a4SJacob Faibussowitsch break; 116a7e14dcfSSatish Balay 117a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 118a7e14dcfSSatish Balay /* Get Reduced Hessian 119a7e14dcfSSatish Balay Msub[i,j] = M[i,j] if i,j in Free_Local or i==j 120a7e14dcfSSatish Balay Msub[i,j] = 0 if i!=j and i or j not in Free_Local 121a7e14dcfSSatish Balay */ 122d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)M); 1239566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-overwrite_hessian", "modify the existing hessian matrix when computing submatrices", "TaoSubsetType", flg, &flg, NULL)); 124d0609cedSBarry Smith PetscOptionsEnd(); 1258afaa268SBarry Smith if (flg) { 1269566063dSJacob Faibussowitsch PetscCall(MatDuplicate(M, MAT_COPY_VALUES, Msub)); 127a7e14dcfSSatish Balay } else { 128a7e14dcfSSatish Balay /* Act on hessian directly (default) */ 1299566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)M)); 130a7e14dcfSSatish Balay *Msub = M; 131a7e14dcfSSatish Balay } 132a7e14dcfSSatish Balay /* Save the diagonal to temporary vector */ 1339566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(*Msub, v1)); 134a7e14dcfSSatish Balay 135a7e14dcfSSatish Balay /* Zero out rows and columns */ 1369566063dSJacob Faibussowitsch PetscCall(ISComplementVec(is, v1, &iscomp)); 137a7e14dcfSSatish Balay 138a7e14dcfSSatish Balay /* Use v1 instead of 0 here because of PETSc bug */ 1399566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsIS(*Msub, iscomp, 1.0, v1, v1)); 140a7e14dcfSSatish Balay 1419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscomp)); 142a7e14dcfSSatish Balay break; 143a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 1449566063dSJacob Faibussowitsch PetscCall(ISComplementVec(is, v1, &iscomp)); 1459566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrixFree(M, iscomp, iscomp, Msub)); 1469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscomp)); 147a7e14dcfSSatish Balay break; 148a7e14dcfSSatish Balay } 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150a7e14dcfSSatish Balay } 1512f75a4aaSAlp Dener 1522f75a4aaSAlp Dener /*@C 1532f75a4aaSAlp Dener TaoEstimateActiveBounds - Generates index sets for variables at the lower and upper 1542f75a4aaSAlp Dener bounds, as well as fixed variables where lower and upper bounds equal each other. 1552f75a4aaSAlp Dener 1562f75a4aaSAlp Dener Input Parameters: 1572f75a4aaSAlp Dener + X - solution vector 1582f75a4aaSAlp Dener . XL - lower bound vector 1592f75a4aaSAlp Dener . XU - upper bound vector 1602f75a4aaSAlp Dener . G - unprojected gradient 1610a4511e9SAlp Dener . S - step direction with which the active bounds will be estimated 162b6a6cedcSAlp Dener . W - work vector of type and size of X 1630a4511e9SAlp Dener - steplen - the step length at which the active bounds will be estimated (needs to be conservative) 1642f75a4aaSAlp Dener 1652f75a4aaSAlp Dener Output Parameters: 16676be6f4fSStefano Zampini + bound_tol - tolerance for the bound estimation 1672f75a4aaSAlp Dener . active_lower - index set for active variables at the lower bound 1682f75a4aaSAlp Dener . active_upper - index set for active variables at the upper bound 1692f75a4aaSAlp Dener . active_fixed - index set for fixed variables 1702f75a4aaSAlp Dener . active - index set for all active variables 171b6a6cedcSAlp Dener - inactive - complementary index set for inactive variables 1720a4511e9SAlp Dener 1730a4511e9SAlp Dener Notes: 1740a4511e9SAlp Dener This estimation is based on Bertsekas' method, with a built in diagonal scaling value of 1.0e-3. 1750a4511e9SAlp Dener 176b6a6cedcSAlp Dener Level: developer 1772f75a4aaSAlp Dener @*/ 178d71ae5a4SJacob Faibussowitsch 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) 179d71ae5a4SJacob Faibussowitsch { 1802f75a4aaSAlp Dener PetscReal wnorm; 18189da521bSAlp Dener PetscReal zero = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0); 182c4b75bccSAlp Dener PetscInt i, n_isl = 0, n_isu = 0, n_isf = 0, n_isa = 0, n_isi = 0; 183c4b75bccSAlp Dener PetscInt N_isl, N_isu, N_isf, N_isa, N_isi; 184c4b75bccSAlp Dener PetscInt n, low, high, nDiff; 185c4b75bccSAlp Dener PetscInt *isl = NULL, *isu = NULL, *isf = NULL, *isa = NULL, *isi = NULL; 1862f75a4aaSAlp Dener const PetscScalar *xl, *xu, *x, *g; 1876519dc0cSAlp Dener MPI_Comm comm = PetscObjectComm((PetscObject)X); 1882f75a4aaSAlp Dener 1892f75a4aaSAlp Dener PetscFunctionBegin; 190c4b75bccSAlp Dener PetscValidHeaderSpecific(X, VEC_CLASSID, 1); 19176be6f4fSStefano Zampini if (XL) PetscValidHeaderSpecific(XL, VEC_CLASSID, 2); 19276be6f4fSStefano Zampini if (XU) PetscValidHeaderSpecific(XU, VEC_CLASSID, 3); 193c4b75bccSAlp Dener PetscValidHeaderSpecific(G, VEC_CLASSID, 4); 194c4b75bccSAlp Dener PetscValidHeaderSpecific(S, VEC_CLASSID, 5); 195c4b75bccSAlp Dener PetscValidHeaderSpecific(W, VEC_CLASSID, 6); 196c4b75bccSAlp Dener 19776be6f4fSStefano Zampini if (XL) PetscCheckSameType(X, 1, XL, 2); 19876be6f4fSStefano Zampini if (XU) PetscCheckSameType(X, 1, XU, 3); 199c4b75bccSAlp Dener PetscCheckSameType(X, 1, G, 4); 200c4b75bccSAlp Dener PetscCheckSameType(X, 1, S, 5); 201c4b75bccSAlp Dener PetscCheckSameType(X, 1, W, 6); 20276be6f4fSStefano Zampini if (XL) PetscCheckSameComm(X, 1, XL, 2); 20376be6f4fSStefano Zampini if (XU) PetscCheckSameComm(X, 1, XU, 3); 204c4b75bccSAlp Dener PetscCheckSameComm(X, 1, G, 4); 205c4b75bccSAlp Dener PetscCheckSameComm(X, 1, S, 5); 206c4b75bccSAlp Dener PetscCheckSameComm(X, 1, W, 6); 20776be6f4fSStefano Zampini if (XL) VecCheckSameSize(X, 1, XL, 2); 20876be6f4fSStefano Zampini if (XU) VecCheckSameSize(X, 1, XU, 3); 209c4b75bccSAlp Dener VecCheckSameSize(X, 1, G, 4); 210c4b75bccSAlp Dener VecCheckSameSize(X, 1, S, 5); 211c4b75bccSAlp Dener VecCheckSameSize(X, 1, W, 6); 212c4b75bccSAlp Dener 2132f75a4aaSAlp Dener /* Update the tolerance for bound detection (this is based on Bertsekas' method) */ 2149566063dSJacob Faibussowitsch PetscCall(VecCopy(X, W)); 2159566063dSJacob Faibussowitsch PetscCall(VecAXPBY(W, steplen, 1.0, S)); 2169566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(W, XL, XU, 0.0, &nDiff, W)); 2179566063dSJacob Faibussowitsch PetscCall(VecAXPBY(W, 1.0, -1.0, X)); 2189566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 2192f75a4aaSAlp Dener *bound_tol = PetscMin(*bound_tol, wnorm); 2202f75a4aaSAlp Dener 22176be6f4fSStefano Zampini /* Clear all index sets */ 22276be6f4fSStefano Zampini PetscCall(ISDestroy(active_lower)); 22376be6f4fSStefano Zampini PetscCall(ISDestroy(active_upper)); 22476be6f4fSStefano Zampini PetscCall(ISDestroy(active_fixed)); 22576be6f4fSStefano Zampini PetscCall(ISDestroy(active)); 22676be6f4fSStefano Zampini PetscCall(ISDestroy(inactive)); 22776be6f4fSStefano Zampini 2289566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &low, &high)); 2299566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 23076be6f4fSStefano Zampini if (!XL && !XU) { 23176be6f4fSStefano Zampini PetscCall(ISCreateStride(comm, n, low, 1, inactive)); 2323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23376be6f4fSStefano Zampini } 2342f75a4aaSAlp Dener if (n > 0) { 2359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XL, &xl)); 2379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XU, &xu)); 2389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(G, &g)); 2392f75a4aaSAlp Dener 2402f75a4aaSAlp Dener /* Loop over variables and categorize the indexes */ 2419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isl)); 2429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isu)); 2439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isf)); 2449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isa)); 2459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isi)); 246c4b75bccSAlp Dener for (i = 0; i < n; ++i) { 24766b4d56fSTodd Munson if (xl[i] == xu[i]) { 248c4b75bccSAlp Dener /* Fixed variables */ 2499371c9d4SSatish Balay isf[n_isf] = low + i; 2509371c9d4SSatish Balay ++n_isf; 2519371c9d4SSatish Balay isa[n_isa] = low + i; 2529371c9d4SSatish Balay ++n_isa; 25376be6f4fSStefano Zampini } else if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + *bound_tol && g[i] > zero) { 254c4b75bccSAlp Dener /* Lower bounded variables */ 2559371c9d4SSatish Balay isl[n_isl] = low + i; 2569371c9d4SSatish Balay ++n_isl; 2579371c9d4SSatish Balay isa[n_isa] = low + i; 2589371c9d4SSatish Balay ++n_isa; 25976be6f4fSStefano Zampini } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - *bound_tol && g[i] < zero) { 260c4b75bccSAlp Dener /* Upper bounded variables */ 2619371c9d4SSatish Balay isu[n_isu] = low + i; 2629371c9d4SSatish Balay ++n_isu; 2639371c9d4SSatish Balay isa[n_isa] = low + i; 2649371c9d4SSatish Balay ++n_isa; 265c4b75bccSAlp Dener } else { 266c4b75bccSAlp Dener /* Inactive variables */ 2679371c9d4SSatish Balay isi[n_isi] = low + i; 2689371c9d4SSatish Balay ++n_isi; 2692f75a4aaSAlp Dener } 2702f75a4aaSAlp Dener } 2712f75a4aaSAlp Dener 2729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XL, &xl)); 2749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XU, &xu)); 2759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(G, &g)); 2762f75a4aaSAlp Dener } 2772f75a4aaSAlp Dener 278c4b75bccSAlp Dener /* Collect global sizes */ 2791c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&n_isl, &N_isl, 1, MPIU_INT, MPI_SUM, comm)); 2801c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&n_isu, &N_isu, 1, MPIU_INT, MPI_SUM, comm)); 2811c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&n_isf, &N_isf, 1, MPIU_INT, MPI_SUM, comm)); 2821c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&n_isa, &N_isa, 1, MPIU_INT, MPI_SUM, comm)); 2831c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&n_isi, &N_isi, 1, MPIU_INT, MPI_SUM, comm)); 284c4b75bccSAlp Dener 285c4b75bccSAlp Dener /* Create index set for lower bounded variables */ 286c4b75bccSAlp Dener if (N_isl > 0) { 2879566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_isl, isl, PETSC_OWN_POINTER, active_lower)); 2887e9273c4SAlp Dener } else { 2899566063dSJacob Faibussowitsch PetscCall(PetscFree(isl)); 290c4b75bccSAlp Dener } 291c4b75bccSAlp Dener /* Create index set for upper bounded variables */ 292c4b75bccSAlp Dener if (N_isu > 0) { 2939566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_isu, isu, PETSC_OWN_POINTER, active_upper)); 2947e9273c4SAlp Dener } else { 2959566063dSJacob Faibussowitsch PetscCall(PetscFree(isu)); 296c4b75bccSAlp Dener } 297c4b75bccSAlp Dener /* Create index set for fixed variables */ 298c4b75bccSAlp Dener if (N_isf > 0) { 2999566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_isf, isf, PETSC_OWN_POINTER, active_fixed)); 3007e9273c4SAlp Dener } else { 3019566063dSJacob Faibussowitsch PetscCall(PetscFree(isf)); 302c4b75bccSAlp Dener } 303c4b75bccSAlp Dener /* Create index set for all actively bounded variables */ 304c4b75bccSAlp Dener if (N_isa > 0) { 3059566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_isa, isa, PETSC_OWN_POINTER, active)); 3067e9273c4SAlp Dener } else { 3079566063dSJacob Faibussowitsch PetscCall(PetscFree(isa)); 308c4b75bccSAlp Dener } 309c4b75bccSAlp Dener /* Create index set for all inactive variables */ 310c4b75bccSAlp Dener if (N_isi > 0) { 3119566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_isi, isi, PETSC_OWN_POINTER, inactive)); 3127e9273c4SAlp Dener } else { 3139566063dSJacob Faibussowitsch PetscCall(PetscFree(isi)); 314c4b75bccSAlp Dener } 3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3162f75a4aaSAlp Dener } 3172f75a4aaSAlp Dener 3182f75a4aaSAlp Dener /*@C 3192f75a4aaSAlp Dener TaoBoundStep - Ensures the correct zero or adjusted step direction 3202f75a4aaSAlp Dener values for active variables. 3212f75a4aaSAlp Dener 3222f75a4aaSAlp Dener Input Parameters: 3232f75a4aaSAlp Dener + X - solution vector 3242f75a4aaSAlp Dener . XL - lower bound vector 3252f75a4aaSAlp Dener . XU - upper bound vector 3262f75a4aaSAlp Dener . active_lower - index set for lower bounded active variables 3272f75a4aaSAlp Dener . active_upper - index set for lower bounded active variables 328b6a6cedcSAlp Dener . active_fixed - index set for fixed active variables 329b6a6cedcSAlp Dener - scale - amplification factor for the step that needs to be taken on actively bounded variables 3302f75a4aaSAlp Dener 33197bb3fdcSJose E. Roman Output Parameter: 3322f75a4aaSAlp Dener . S - step direction to be modified 333b6a6cedcSAlp Dener 334b6a6cedcSAlp Dener Level: developer 3352f75a4aaSAlp Dener @*/ 336d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, PetscReal scale, Vec S) 337d71ae5a4SJacob Faibussowitsch { 3382f75a4aaSAlp Dener Vec step_lower, step_upper, step_fixed; 3392f75a4aaSAlp Dener Vec x_lower, x_upper; 3402f75a4aaSAlp Dener Vec bound_lower, bound_upper; 3412f75a4aaSAlp Dener 3422f75a4aaSAlp Dener PetscFunctionBegin; 3432f75a4aaSAlp Dener /* Adjust step for variables at the estimated lower bound */ 3442f75a4aaSAlp Dener if (active_lower) { 3459566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(S, active_lower, &step_lower)); 3469566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, active_lower, &x_lower)); 3479566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(XL, active_lower, &bound_lower)); 3489566063dSJacob Faibussowitsch PetscCall(VecCopy(bound_lower, step_lower)); 3499566063dSJacob Faibussowitsch PetscCall(VecAXPY(step_lower, -1.0, x_lower)); 3509566063dSJacob Faibussowitsch PetscCall(VecScale(step_lower, scale)); 3519566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(S, active_lower, &step_lower)); 3529566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, active_lower, &x_lower)); 3539566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(XL, active_lower, &bound_lower)); 3542f75a4aaSAlp Dener } 3552f75a4aaSAlp Dener 3562f75a4aaSAlp Dener /* Adjust step for the variables at the estimated upper bound */ 3572f75a4aaSAlp Dener if (active_upper) { 3589566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(S, active_upper, &step_upper)); 3599566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, active_upper, &x_upper)); 3609566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(XU, active_upper, &bound_upper)); 3619566063dSJacob Faibussowitsch PetscCall(VecCopy(bound_upper, step_upper)); 3629566063dSJacob Faibussowitsch PetscCall(VecAXPY(step_upper, -1.0, x_upper)); 3639566063dSJacob Faibussowitsch PetscCall(VecScale(step_upper, scale)); 3649566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(S, active_upper, &step_upper)); 3659566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, active_upper, &x_upper)); 3669566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(XU, active_upper, &bound_upper)); 3672f75a4aaSAlp Dener } 3682f75a4aaSAlp Dener 3692f75a4aaSAlp Dener /* Zero out step for fixed variables */ 3702f75a4aaSAlp Dener if (active_fixed) { 3719566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(S, active_fixed, &step_fixed)); 3729566063dSJacob Faibussowitsch PetscCall(VecSet(step_fixed, 0.0)); 3739566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(S, active_fixed, &step_fixed)); 3742f75a4aaSAlp Dener } 3753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3762f75a4aaSAlp Dener } 377c4b75bccSAlp Dener 378c4b75bccSAlp Dener /*@C 379b6a6cedcSAlp Dener TaoBoundSolution - Ensures that the solution vector is snapped into the bounds within a given tolerance. 380c4b75bccSAlp Dener 381c3339decSBarry Smith Collective 382f0fc11ceSJed Brown 383c4b75bccSAlp Dener Input Parameters: 384b6a6cedcSAlp Dener + X - solution vector 385b6a6cedcSAlp Dener . XL - lower bound vector 386c4b75bccSAlp Dener . XU - upper bound vector 387f0fc11ceSJed Brown - bound_tol - absolute tolerance in enforcing the bound 388c4b75bccSAlp Dener 389c4b75bccSAlp Dener Output Parameters: 390f0fc11ceSJed Brown + nDiff - total number of vector entries that have been bounded 391f0fc11ceSJed Brown - Xout - modified solution vector satisfying bounds to bound_tol 392b6a6cedcSAlp Dener 393b6a6cedcSAlp Dener Level: developer 394f0fc11ceSJed Brown 395db781477SPatrick Sanan .seealso: `TAOBNCG`, `TAOBNTL`, `TAOBNTR` 396c4b75bccSAlp Dener @*/ 397d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBoundSolution(Vec X, Vec XL, Vec XU, PetscReal bound_tol, PetscInt *nDiff, Vec Xout) 398d71ae5a4SJacob Faibussowitsch { 399c4b75bccSAlp Dener PetscInt i, n, low, high, nDiff_loc = 0; 4003b063059SAlp Dener PetscScalar *xout; 4013b063059SAlp Dener const PetscScalar *x, *xl, *xu; 402c4b75bccSAlp Dener 403c4b75bccSAlp Dener PetscFunctionBegin; 404c4b75bccSAlp Dener PetscValidHeaderSpecific(X, VEC_CLASSID, 1); 40576be6f4fSStefano Zampini if (XL) PetscValidHeaderSpecific(XL, VEC_CLASSID, 2); 40676be6f4fSStefano Zampini if (XU) PetscValidHeaderSpecific(XU, VEC_CLASSID, 3); 407064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(Xout, VEC_CLASSID, 6); 40876be6f4fSStefano Zampini if (!XL && !XU) { 40976be6f4fSStefano Zampini PetscCall(VecCopy(X, Xout)); 41076be6f4fSStefano Zampini *nDiff = 0.0; 4113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41276be6f4fSStefano Zampini } 413c4b75bccSAlp Dener PetscCheckSameType(X, 1, XL, 2); 414c4b75bccSAlp Dener PetscCheckSameType(X, 1, XU, 3); 415064a246eSJacob Faibussowitsch PetscCheckSameType(X, 1, Xout, 6); 416c4b75bccSAlp Dener PetscCheckSameComm(X, 1, XL, 2); 417c4b75bccSAlp Dener PetscCheckSameComm(X, 1, XU, 3); 418064a246eSJacob Faibussowitsch PetscCheckSameComm(X, 1, Xout, 6); 419c4b75bccSAlp Dener VecCheckSameSize(X, 1, XL, 2); 420c4b75bccSAlp Dener VecCheckSameSize(X, 1, XU, 3); 4213b063059SAlp Dener VecCheckSameSize(X, 1, Xout, 4); 422c4b75bccSAlp Dener 4239566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &low, &high)); 4249566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 425c4b75bccSAlp Dener if (n > 0) { 4269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 4279566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XL, &xl)); 4289566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XU, &xu)); 4299566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xout, &xout)); 430c4b75bccSAlp Dener 431c4b75bccSAlp Dener for (i = 0; i < n; ++i) { 43276be6f4fSStefano Zampini if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + bound_tol) { 4339371c9d4SSatish Balay xout[i] = xl[i]; 4349371c9d4SSatish Balay ++nDiff_loc; 43576be6f4fSStefano Zampini } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - bound_tol) { 4369371c9d4SSatish Balay xout[i] = xu[i]; 4379371c9d4SSatish Balay ++nDiff_loc; 438c4b75bccSAlp Dener } 439c4b75bccSAlp Dener } 440c4b75bccSAlp Dener 4419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 4429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XL, &xl)); 4439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XU, &xu)); 4449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xout, &xout)); 445c4b75bccSAlp Dener } 4461c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&nDiff_loc, nDiff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)X))); 4473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 448c4b75bccSAlp Dener } 449