1ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/taoimpl.h> 3aaa7dc30SBarry Smith #include <../src/tao/matrix/submatfree.h> 4a7e14dcfSSatish Balay 5b98f30f2SJason Sarich /*@C 6b98f30f2SJason Sarich TaoVecGetSubVec - Gets a subvector using the IS 7a7e14dcfSSatish Balay 8a7e14dcfSSatish Balay Input Parameters: 9a7e14dcfSSatish Balay + vfull - the full matrix 10a7e14dcfSSatish Balay . is - the index set for the subvector 11a7e14dcfSSatish Balay . reduced_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK, TAO_SUBSET_MATRIXFREE) 12a7e14dcfSSatish Balay - maskvalue - the value to set the unused vector elements to (for TAO_SUBSET_MASK or TAO_SUBSET_MATRIXFREE) 13a7e14dcfSSatish Balay 14a7e14dcfSSatish Balay Output Parameters: 15a7e14dcfSSatish Balay . vreduced - the subvector 16a7e14dcfSSatish Balay 171eb8069cSJason Sarich Notes: 18a7e14dcfSSatish Balay maskvalue should usually be 0.0, unless a pointwise divide will be used. 191eb8069cSJason Sarich 20a7e14dcfSSatish Balay @*/ 213a831ad5SBarry Smith PetscErrorCode TaoVecGetSubVec(Vec vfull, IS is, TaoSubsetType reduced_type, PetscReal maskvalue, Vec *vreduced) 22a7e14dcfSSatish Balay { 23a7e14dcfSSatish Balay PetscErrorCode ierr; 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 37a7e14dcfSSatish Balay ierr = VecGetSize(vfull, &nfull);CHKERRQ(ierr); 38a7e14dcfSSatish Balay ierr = ISGetSize(is, &nreduced);CHKERRQ(ierr); 39a7e14dcfSSatish Balay 40a7e14dcfSSatish Balay if (nreduced == nfull) { 41a7e14dcfSSatish Balay ierr = VecDestroy(vreduced);CHKERRQ(ierr); 42a7e14dcfSSatish Balay ierr = VecDuplicate(vfull,vreduced);CHKERRQ(ierr); 43a7e14dcfSSatish Balay ierr = VecCopy(vfull,*vreduced);CHKERRQ(ierr); 44a7e14dcfSSatish Balay } else { 45a7e14dcfSSatish Balay switch (reduced_type) { 46a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 47a7e14dcfSSatish Balay ierr = VecGetType(vfull,&vtype);CHKERRQ(ierr); 48a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(vfull,&flow,&fhigh);CHKERRQ(ierr); 49a7e14dcfSSatish Balay ierr = ISGetLocalSize(is,&nreduced_local);CHKERRQ(ierr); 50a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)vfull,&comm);CHKERRQ(ierr); 51a7e14dcfSSatish Balay if (*vreduced) { 52a7e14dcfSSatish Balay ierr = VecDestroy(vreduced);CHKERRQ(ierr); 53a7e14dcfSSatish Balay } 54a7e14dcfSSatish Balay ierr = VecCreate(comm,vreduced);CHKERRQ(ierr); 55a7e14dcfSSatish Balay ierr = VecSetType(*vreduced,vtype);CHKERRQ(ierr); 56a7e14dcfSSatish Balay 57a7e14dcfSSatish Balay ierr = VecSetSizes(*vreduced,nreduced_local,nreduced);CHKERRQ(ierr); 58a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(*vreduced,&rlow,&rhigh);CHKERRQ(ierr); 59a7e14dcfSSatish Balay ierr = ISCreateStride(comm,nreduced_local,rlow,1,&ident);CHKERRQ(ierr); 60a7e14dcfSSatish Balay ierr = VecScatterCreate(vfull,is,*vreduced,ident,&scatter);CHKERRQ(ierr); 61a7e14dcfSSatish Balay ierr = VecScatterBegin(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 62a7e14dcfSSatish Balay ierr = VecScatterEnd(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 63a7e14dcfSSatish Balay ierr = VecScatterDestroy(&scatter);CHKERRQ(ierr); 64a7e14dcfSSatish Balay ierr = ISDestroy(&ident);CHKERRQ(ierr); 65a7e14dcfSSatish Balay break; 66a7e14dcfSSatish Balay 67a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 68a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 69a7e14dcfSSatish Balay /* vr[i] = vf[i] if i in is 70a7e14dcfSSatish Balay vr[i] = 0 otherwise */ 716c4ed002SBarry Smith if (!*vreduced) { 72a7e14dcfSSatish Balay ierr = VecDuplicate(vfull,vreduced);CHKERRQ(ierr); 73a7e14dcfSSatish Balay } 74a7e14dcfSSatish Balay 75a7e14dcfSSatish Balay ierr = VecSet(*vreduced,maskvalue);CHKERRQ(ierr); 76a7e14dcfSSatish Balay ierr = ISGetLocalSize(is,&nlocal);CHKERRQ(ierr); 77a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(vfull,&flow,&fhigh);CHKERRQ(ierr); 78a7e14dcfSSatish Balay ierr = VecGetArray(vfull,&fv);CHKERRQ(ierr); 79a7e14dcfSSatish Balay ierr = VecGetArray(*vreduced,&rv);CHKERRQ(ierr); 80a7e14dcfSSatish Balay ierr = ISGetIndices(is,&s);CHKERRQ(ierr); 8153506e15SBarry Smith if (nlocal > (fhigh-flow)) SETERRQ2(PETSC_COMM_WORLD,1,"IS local size %d > Vec local size %d",nlocal,fhigh-flow); 82a7e14dcfSSatish Balay for (i=0;i<nlocal;i++) { 83a7e14dcfSSatish Balay rv[s[i]-flow] = fv[s[i]-flow]; 84a7e14dcfSSatish Balay } 85a7e14dcfSSatish Balay ierr = ISRestoreIndices(is,&s);CHKERRQ(ierr); 86a7e14dcfSSatish Balay ierr = VecRestoreArray(vfull,&fv);CHKERRQ(ierr); 87a7e14dcfSSatish Balay ierr = VecRestoreArray(*vreduced,&rv);CHKERRQ(ierr); 88a7e14dcfSSatish Balay break; 89a7e14dcfSSatish Balay } 90a7e14dcfSSatish Balay } 91a7e14dcfSSatish Balay PetscFunctionReturn(0); 92a7e14dcfSSatish Balay } 93a7e14dcfSSatish Balay 94b98f30f2SJason Sarich /*@C 95b98f30f2SJason Sarich TaoMatGetSubMat - Gets a submatrix using the IS 96a7e14dcfSSatish Balay 97a7e14dcfSSatish Balay Input Parameters: 98a7e14dcfSSatish Balay + M - the full matrix (n x n) 99a7e14dcfSSatish Balay . is - the index set for the submatrix (both row and column index sets need to be the same) 100a7e14dcfSSatish Balay . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option 101a7e14dcfSSatish Balay - subset_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK, 102a7e14dcfSSatish Balay TAO_SUBSET_MATRIXFREE) 103a7e14dcfSSatish Balay 104a7e14dcfSSatish Balay Output Parameters: 105a7e14dcfSSatish Balay . Msub - the submatrix 106a7e14dcfSSatish Balay @*/ 107b98f30f2SJason Sarich PetscErrorCode TaoMatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub) 108a7e14dcfSSatish Balay { 109a7e14dcfSSatish Balay PetscErrorCode ierr; 110a7e14dcfSSatish Balay IS iscomp; 111*62675beeSAlp Dener PetscBool flg = PETSC_TRUE; 11253506e15SBarry Smith 113a7e14dcfSSatish Balay PetscFunctionBegin; 114a7e14dcfSSatish Balay PetscValidHeaderSpecific(M,MAT_CLASSID,1); 115a7e14dcfSSatish Balay PetscValidHeaderSpecific(is,IS_CLASSID,2); 116a7e14dcfSSatish Balay ierr = MatDestroy(Msub);CHKERRQ(ierr); 117a7e14dcfSSatish Balay switch (subset_type) { 118a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 1197dae84e0SHong Zhang ierr = MatCreateSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub);CHKERRQ(ierr); 120a7e14dcfSSatish Balay break; 121a7e14dcfSSatish Balay 122a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 123a7e14dcfSSatish Balay /* Get Reduced Hessian 124a7e14dcfSSatish Balay Msub[i,j] = M[i,j] if i,j in Free_Local or i==j 125a7e14dcfSSatish Balay Msub[i,j] = 0 if i!=j and i or j not in Free_Local 126a7e14dcfSSatish Balay */ 1271a1499c8SBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)M),NULL,NULL,NULL);CHKERRQ(ierr); 128*62675beeSAlp Dener ierr = PetscOptionsBool("-overwrite_hessian","modify the existing hessian matrix when computing submatrices","TaoSubsetType",flg,&flg,NULL);CHKERRQ(ierr); 129302440fdSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 1308afaa268SBarry Smith if (flg) { 131a7e14dcfSSatish Balay ierr = MatDuplicate(M, MAT_COPY_VALUES, Msub);CHKERRQ(ierr); 132a7e14dcfSSatish Balay } else { 133a7e14dcfSSatish Balay /* Act on hessian directly (default) */ 134a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)M);CHKERRQ(ierr); 135a7e14dcfSSatish Balay *Msub = M; 136a7e14dcfSSatish Balay } 137a7e14dcfSSatish Balay /* Save the diagonal to temporary vector */ 138a7e14dcfSSatish Balay ierr = MatGetDiagonal(*Msub,v1);CHKERRQ(ierr); 139a7e14dcfSSatish Balay 140a7e14dcfSSatish Balay /* Zero out rows and columns */ 1414473680cSBarry Smith ierr = ISComplementVec(is,v1,&iscomp);CHKERRQ(ierr); 142a7e14dcfSSatish Balay 143a7e14dcfSSatish Balay /* Use v1 instead of 0 here because of PETSc bug */ 144a7e14dcfSSatish Balay ierr = MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1);CHKERRQ(ierr); 145a7e14dcfSSatish Balay 146a7e14dcfSSatish Balay ierr = ISDestroy(&iscomp);CHKERRQ(ierr); 147a7e14dcfSSatish Balay break; 148a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 1494473680cSBarry Smith ierr = ISComplementVec(is,v1,&iscomp);CHKERRQ(ierr); 150a7e14dcfSSatish Balay ierr = MatCreateSubMatrixFree(M,iscomp,iscomp,Msub);CHKERRQ(ierr); 151a7e14dcfSSatish Balay ierr = ISDestroy(&iscomp);CHKERRQ(ierr); 152a7e14dcfSSatish Balay break; 153a7e14dcfSSatish Balay } 154a7e14dcfSSatish Balay PetscFunctionReturn(0); 155a7e14dcfSSatish Balay } 1562f75a4aaSAlp Dener 1572f75a4aaSAlp Dener /*@C 1582f75a4aaSAlp Dener TaoEstimateActiveBounds - Generates index sets for variables at the lower and upper 1592f75a4aaSAlp Dener bounds, as well as fixed variables where lower and upper bounds equal each other. 1602f75a4aaSAlp Dener 1612f75a4aaSAlp Dener Input Parameters: 1622f75a4aaSAlp Dener + X - solution vector 1632f75a4aaSAlp Dener . XL - lower bound vector 1642f75a4aaSAlp Dener . XU - upper bound vector 1652f75a4aaSAlp Dener . G - unprojected gradient 1660a4511e9SAlp Dener . S - step direction with which the active bounds will be estimated 1670a4511e9SAlp Dener - steplen - the step length at which the active bounds will be estimated (needs to be conservative) 1682f75a4aaSAlp Dener 1692f75a4aaSAlp Dener Output Parameters: 1702f75a4aaSAlp Dener . bound_tol - tolerance for for the bound estimation 1712f75a4aaSAlp Dener . active_lower - index set for active variables at the lower bound 1722f75a4aaSAlp Dener . active_upper - index set for active variables at the upper bound 1732f75a4aaSAlp Dener . active_fixed - index set for fixed variables 1742f75a4aaSAlp Dener . active - index set for all active variables 1752f75a4aaSAlp Dener . inactive - complementary index set for inactive variables 1760a4511e9SAlp Dener 1770a4511e9SAlp Dener Notes: 1780a4511e9SAlp Dener This estimation is based on Bertsekas' method, with a built in diagonal scaling value of 1.0e-3. 1790a4511e9SAlp Dener 1802f75a4aaSAlp Dener @*/ 1810a4511e9SAlp Dener PetscErrorCode TaoEstimateActiveBounds(Vec X, Vec XL, Vec XU, Vec G, Vec S, PetscReal steplen, PetscReal *bound_tol, 1822f75a4aaSAlp Dener IS *active_lower, IS *active_upper, IS *active_fixed, IS *active, IS *inactive) 1832f75a4aaSAlp Dener { 1842f75a4aaSAlp Dener PetscErrorCode ierr; 1852f75a4aaSAlp Dener 1862f75a4aaSAlp Dener Vec W; 1872f75a4aaSAlp Dener PetscReal wnorm; 1882f75a4aaSAlp Dener PetscInt i, n_isl=0, n_isu=0, n_isf=0; 1892f75a4aaSAlp Dener PetscInt n, low, high; 1902f75a4aaSAlp Dener PetscInt *isl=NULL, *isu=NULL, *isf=NULL; 1912f75a4aaSAlp Dener const PetscScalar *xl, *xu, *x, *g; 1922f75a4aaSAlp Dener 1932f75a4aaSAlp Dener PetscFunctionBegin; 1942f75a4aaSAlp Dener /* Update the tolerance for bound detection (this is based on Bertsekas' method) */ 1952f75a4aaSAlp Dener ierr = VecDuplicate(S, &W);CHKERRQ(ierr); 1962f75a4aaSAlp Dener ierr = VecCopy(S, W);CHKERRQ(ierr); 1970a4511e9SAlp Dener ierr = VecAXPBY(W, 1.0, steplen, X);CHKERRQ(ierr); 1982f75a4aaSAlp Dener ierr = VecMedian(XL, W, XU, W);CHKERRQ(ierr); 1992f75a4aaSAlp Dener ierr = VecAXPBY(W, 1.0, -1.0, X);CHKERRQ(ierr); 2002f75a4aaSAlp Dener ierr = VecNorm(W, NORM_2, &wnorm);CHKERRQ(ierr); 2012f75a4aaSAlp Dener *bound_tol = PetscMin(*bound_tol, wnorm); 2022f75a4aaSAlp Dener 2032f75a4aaSAlp Dener ierr = VecGetOwnershipRange(X, &low, &high);CHKERRQ(ierr); 2042f75a4aaSAlp Dener ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr); 2052f75a4aaSAlp Dener if (n>0){ 2062f75a4aaSAlp Dener ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 2072f75a4aaSAlp Dener ierr = VecGetArrayRead(XL, &xl);CHKERRQ(ierr); 2082f75a4aaSAlp Dener ierr = VecGetArrayRead(XU, &xu);CHKERRQ(ierr); 2092f75a4aaSAlp Dener ierr = VecGetArrayRead(G, &g);CHKERRQ(ierr); 2102f75a4aaSAlp Dener 2112f75a4aaSAlp Dener /* Loop over variables and categorize the indexes */ 2122f75a4aaSAlp Dener ierr = PetscMalloc1(n, &isl);CHKERRQ(ierr); 2132f75a4aaSAlp Dener ierr = PetscMalloc1(n, &isu);CHKERRQ(ierr); 2142f75a4aaSAlp Dener ierr = PetscMalloc1(n, &isf);CHKERRQ(ierr); 2152f75a4aaSAlp Dener for (i=0; i<n; i++) { 2162f75a4aaSAlp Dener if (xl[i] == xu[i]) { 2172f75a4aaSAlp Dener /* Fixed variables here */ 2182f75a4aaSAlp Dener isf[n_isf]=low+i; ++n_isf; 2192f75a4aaSAlp Dener } else if ((x[i] <= xl[i] + *bound_tol) && (g[i] > 0.0)) { 2202f75a4aaSAlp Dener /* Lower bounded variables here */ 2212f75a4aaSAlp Dener isl[n_isl]=low+i; ++n_isl; 2222f75a4aaSAlp Dener } else if ((x[i] >= xu[i] - *bound_tol) && (g[i] < 0.0)) { 2232f75a4aaSAlp Dener /* Upper bounded variables here */ 2242f75a4aaSAlp Dener isu[n_isu]=low+i; ++n_isu; 2252f75a4aaSAlp Dener } 2262f75a4aaSAlp Dener } 2272f75a4aaSAlp Dener 2282f75a4aaSAlp Dener ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 2292f75a4aaSAlp Dener ierr = VecRestoreArrayRead(XL, &xl);CHKERRQ(ierr); 2302f75a4aaSAlp Dener ierr = VecRestoreArrayRead(XU, &xu);CHKERRQ(ierr); 2312f75a4aaSAlp Dener ierr = VecRestoreArrayRead(G, &g);CHKERRQ(ierr); 2322f75a4aaSAlp Dener } 2332f75a4aaSAlp Dener 2342f75a4aaSAlp Dener /* Create index set for lower bounded variables */ 2352f75a4aaSAlp Dener ierr = ISDestroy(active_lower);CHKERRQ(ierr); 2362f75a4aaSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isl, isl, PETSC_OWN_POINTER, active_lower);CHKERRQ(ierr); 2372f75a4aaSAlp Dener /* Create index set for upper bounded variables */ 2382f75a4aaSAlp Dener ierr = ISDestroy(active_upper);CHKERRQ(ierr); 2392f75a4aaSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isu, isu, PETSC_OWN_POINTER, active_upper);CHKERRQ(ierr); 2402f75a4aaSAlp Dener /* Create index set for fixed variables */ 2412f75a4aaSAlp Dener ierr = ISDestroy(active_fixed);CHKERRQ(ierr); 2422f75a4aaSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isf, isf, PETSC_OWN_POINTER, active_fixed);CHKERRQ(ierr); 2432f75a4aaSAlp Dener 2442f75a4aaSAlp Dener /* Create the combined active set */ 2452f75a4aaSAlp Dener ierr = ISDestroy(active);CHKERRQ(ierr); 2462f75a4aaSAlp Dener if (*active_lower && *active_upper && *active_fixed) { 2472f75a4aaSAlp Dener /* All three types of active variables exist */ 2482f75a4aaSAlp Dener const IS islist[3] = {*active_lower, *active_upper, *active_fixed}; 2492f75a4aaSAlp Dener ierr = ISConcatenate(PetscObjectComm((PetscObject)X), 3, islist, active);CHKERRQ(ierr); 2502f75a4aaSAlp Dener ierr = ISSort(*active);CHKERRQ(ierr); 2512f75a4aaSAlp Dener } else if (*active_lower && *active_upper) { 2522f75a4aaSAlp Dener /* Only lower and upper bounded active variables exist */ 2532f75a4aaSAlp Dener ierr = ISSum(*active_lower, *active_upper, active);CHKERRQ(ierr); 2542f75a4aaSAlp Dener } else if (*active_lower && *active_fixed) { 2552f75a4aaSAlp Dener /* Only lower bounded and fixed active variables exist */ 2562f75a4aaSAlp Dener ierr = ISSum(*active_lower, *active_fixed, active);CHKERRQ(ierr); 2572f75a4aaSAlp Dener } else if (*active_upper && *active_fixed) { 2582f75a4aaSAlp Dener /* Only upper bounded and fixed active variables exist */ 2592f75a4aaSAlp Dener ierr = ISSum(*active_upper, *active_fixed, active);CHKERRQ(ierr); 2602f75a4aaSAlp Dener } else if (*active_lower) { 2612f75a4aaSAlp Dener /* Only lower bounded active variables exist */ 2622f75a4aaSAlp Dener *active = *active_lower; 2632f75a4aaSAlp Dener } else if (*active_upper) { 2642f75a4aaSAlp Dener /* Only upper bounded active variables exist */ 2652f75a4aaSAlp Dener *active = *active_upper; 2662f75a4aaSAlp Dener } else if (*active_fixed) { 2672f75a4aaSAlp Dener /* Only fixed active variables exist */ 2682f75a4aaSAlp Dener *active = *active_fixed; 2692f75a4aaSAlp Dener } 2702f75a4aaSAlp Dener /* Create the inactive set */ 2712f75a4aaSAlp Dener ierr = ISDestroy(inactive);CHKERRQ(ierr); 2722f75a4aaSAlp Dener if (*active) { ierr = ISComplementVec(*active, X, inactive);CHKERRQ(ierr); } 2732f75a4aaSAlp Dener 2742f75a4aaSAlp Dener /* Clean up and exit */ 2752f75a4aaSAlp Dener ierr = VecDestroy(&W);CHKERRQ(ierr); 2762f75a4aaSAlp Dener PetscFunctionReturn(0); 2772f75a4aaSAlp Dener } 2782f75a4aaSAlp Dener 2792f75a4aaSAlp Dener /*@C 2802f75a4aaSAlp Dener TaoBoundStep - Ensures the correct zero or adjusted step direction 2812f75a4aaSAlp Dener values for active variables. 2822f75a4aaSAlp Dener 2832f75a4aaSAlp Dener Input Parameters: 2842f75a4aaSAlp Dener + X - solution vector 2852f75a4aaSAlp Dener . XL - lower bound vector 2862f75a4aaSAlp Dener . XU - upper bound vector 2872f75a4aaSAlp Dener . active_lower - index set for lower bounded active variables 2882f75a4aaSAlp Dener . active_upper - index set for lower bounded active variables 2892f75a4aaSAlp Dener - active_fixed - index set for fixed active variables 2902f75a4aaSAlp Dener 2912f75a4aaSAlp Dener Output Parameters: 2922f75a4aaSAlp Dener . S - step direction to be modified 2932f75a4aaSAlp Dener @*/ 2942f75a4aaSAlp Dener PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, Vec S) 2952f75a4aaSAlp Dener { 2962f75a4aaSAlp Dener PetscErrorCode ierr; 2972f75a4aaSAlp Dener 2982f75a4aaSAlp Dener Vec step_lower, step_upper, step_fixed; 2992f75a4aaSAlp Dener Vec x_lower, x_upper; 3002f75a4aaSAlp Dener Vec bound_lower, bound_upper; 3012f75a4aaSAlp Dener 3022f75a4aaSAlp Dener PetscFunctionBegin; 3032f75a4aaSAlp Dener /* Adjust step for variables at the estimated lower bound */ 3042f75a4aaSAlp Dener if (active_lower) { 3052f75a4aaSAlp Dener ierr = VecGetSubVector(S, active_lower, &step_lower);CHKERRQ(ierr); 3062f75a4aaSAlp Dener ierr = VecGetSubVector(X, active_lower, &x_lower);CHKERRQ(ierr); 3072f75a4aaSAlp Dener ierr = VecGetSubVector(XL, active_lower, &bound_lower);CHKERRQ(ierr); 3082f75a4aaSAlp Dener ierr = VecCopy(bound_lower, step_lower);CHKERRQ(ierr); 3092f75a4aaSAlp Dener ierr = VecAXPY(step_lower, -1.0, x_lower);CHKERRQ(ierr); 3102f75a4aaSAlp Dener ierr = VecRestoreSubVector(S, active_lower, &step_lower);CHKERRQ(ierr); 3112f75a4aaSAlp Dener ierr = VecRestoreSubVector(X, active_lower, &x_lower);CHKERRQ(ierr); 3122f75a4aaSAlp Dener ierr = VecRestoreSubVector(XL, active_lower, &bound_lower);CHKERRQ(ierr); 3132f75a4aaSAlp Dener } 3142f75a4aaSAlp Dener 3152f75a4aaSAlp Dener /* Adjust step for the variables at the estimated upper bound */ 3162f75a4aaSAlp Dener if (active_upper) { 3172f75a4aaSAlp Dener ierr = VecGetSubVector(S, active_upper, &step_upper);CHKERRQ(ierr); 3182f75a4aaSAlp Dener ierr = VecGetSubVector(X, active_upper, &x_upper);CHKERRQ(ierr); 3192f75a4aaSAlp Dener ierr = VecGetSubVector(XU, active_upper, &bound_upper);CHKERRQ(ierr); 3202f75a4aaSAlp Dener ierr = VecCopy(bound_upper, step_upper);CHKERRQ(ierr); 3212f75a4aaSAlp Dener ierr = VecAXPY(step_upper, -1.0, x_upper);CHKERRQ(ierr); 3222f75a4aaSAlp Dener ierr = VecRestoreSubVector(S, active_upper, &step_upper);CHKERRQ(ierr); 3232f75a4aaSAlp Dener ierr = VecRestoreSubVector(X, active_upper, &x_upper);CHKERRQ(ierr); 3242f75a4aaSAlp Dener ierr = VecRestoreSubVector(XU, active_upper, &bound_upper);CHKERRQ(ierr); 3252f75a4aaSAlp Dener } 3262f75a4aaSAlp Dener 3272f75a4aaSAlp Dener /* Zero out step for fixed variables */ 3282f75a4aaSAlp Dener if (active_fixed) { 3292f75a4aaSAlp Dener ierr = VecGetSubVector(S, active_fixed, &step_fixed);CHKERRQ(ierr); 3302f75a4aaSAlp Dener ierr = VecSet(step_fixed, 0.0);CHKERRQ(ierr); 3312f75a4aaSAlp Dener ierr = VecRestoreSubVector(S, active_fixed, &step_fixed);CHKERRQ(ierr); 3322f75a4aaSAlp Dener } 3332f75a4aaSAlp Dener PetscFunctionReturn(0); 3342f75a4aaSAlp Dener } 335