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; 1118afaa268SBarry Smith PetscBool flg = PETSC_FALSE; 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); 128302440fdSBarry Smith ierr = PetscOptionsBool("-different_submatrix","use separate 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 } 156*2f75a4aaSAlp Dener 157*2f75a4aaSAlp Dener /*@C 158*2f75a4aaSAlp Dener TaoEstimateActiveBounds - Generates index sets for variables at the lower and upper 159*2f75a4aaSAlp Dener bounds, as well as fixed variables where lower and upper bounds equal each other. 160*2f75a4aaSAlp Dener 161*2f75a4aaSAlp Dener Input Parameters: 162*2f75a4aaSAlp Dener + X - solution vector 163*2f75a4aaSAlp Dener . XL - lower bound vector 164*2f75a4aaSAlp Dener . XU - upper bound vector 165*2f75a4aaSAlp Dener . G - unprojected gradient 166*2f75a4aaSAlp Dener - S - step direction with which the active bounds will be estimated 167*2f75a4aaSAlp Dener 168*2f75a4aaSAlp Dener Output Parameters: 169*2f75a4aaSAlp Dener . bound_tol - tolerance for for the bound estimation 170*2f75a4aaSAlp Dener . active_lower - index set for active variables at the lower bound 171*2f75a4aaSAlp Dener . active_upper - index set for active variables at the upper bound 172*2f75a4aaSAlp Dener . active_fixed - index set for fixed variables 173*2f75a4aaSAlp Dener . active - index set for all active variables 174*2f75a4aaSAlp Dener . inactive - complementary index set for inactive variables 175*2f75a4aaSAlp Dener @*/ 176*2f75a4aaSAlp Dener PetscErrorCode TaoEstimateActiveBounds(Vec X, Vec XL, Vec XU, Vec G, Vec S, PetscReal *bound_tol, 177*2f75a4aaSAlp Dener IS *active_lower, IS *active_upper, IS *active_fixed, IS *active, IS *inactive) 178*2f75a4aaSAlp Dener { 179*2f75a4aaSAlp Dener PetscErrorCode ierr; 180*2f75a4aaSAlp Dener 181*2f75a4aaSAlp Dener Vec W; 182*2f75a4aaSAlp Dener PetscReal wnorm; 183*2f75a4aaSAlp Dener PetscInt i, n_isl=0, n_isu=0, n_isf=0; 184*2f75a4aaSAlp Dener PetscInt n, low, high; 185*2f75a4aaSAlp Dener PetscInt *isl=NULL, *isu=NULL, *isf=NULL; 186*2f75a4aaSAlp Dener const PetscScalar *xl, *xu, *x, *g; 187*2f75a4aaSAlp Dener 188*2f75a4aaSAlp Dener PetscFunctionBegin; 189*2f75a4aaSAlp Dener /* Update the tolerance for bound detection (this is based on Bertsekas' method) */ 190*2f75a4aaSAlp Dener ierr = VecDuplicate(S, &W);CHKERRQ(ierr); 191*2f75a4aaSAlp Dener ierr = VecCopy(S, W);CHKERRQ(ierr); 192*2f75a4aaSAlp Dener ierr = VecAXPBY(W, 1.0, 0.001, X);CHKERRQ(ierr); 193*2f75a4aaSAlp Dener ierr = VecMedian(XL, W, XU, W);CHKERRQ(ierr); 194*2f75a4aaSAlp Dener ierr = VecAXPBY(W, 1.0, -1.0, X);CHKERRQ(ierr); 195*2f75a4aaSAlp Dener ierr = VecNorm(W, NORM_2, &wnorm);CHKERRQ(ierr); 196*2f75a4aaSAlp Dener *bound_tol = PetscMin(*bound_tol, wnorm); 197*2f75a4aaSAlp Dener 198*2f75a4aaSAlp Dener ierr = VecGetOwnershipRange(X, &low, &high);CHKERRQ(ierr); 199*2f75a4aaSAlp Dener ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr); 200*2f75a4aaSAlp Dener if (n>0){ 201*2f75a4aaSAlp Dener ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 202*2f75a4aaSAlp Dener ierr = VecGetArrayRead(XL, &xl);CHKERRQ(ierr); 203*2f75a4aaSAlp Dener ierr = VecGetArrayRead(XU, &xu);CHKERRQ(ierr); 204*2f75a4aaSAlp Dener ierr = VecGetArrayRead(G, &g);CHKERRQ(ierr); 205*2f75a4aaSAlp Dener 206*2f75a4aaSAlp Dener /* Loop over variables and categorize the indexes */ 207*2f75a4aaSAlp Dener ierr = PetscMalloc1(n, &isl);CHKERRQ(ierr); 208*2f75a4aaSAlp Dener ierr = PetscMalloc1(n, &isu);CHKERRQ(ierr); 209*2f75a4aaSAlp Dener ierr = PetscMalloc1(n, &isf);CHKERRQ(ierr); 210*2f75a4aaSAlp Dener for (i=0; i<n; i++) { 211*2f75a4aaSAlp Dener if (xl[i] == xu[i]) { 212*2f75a4aaSAlp Dener /* Fixed variables here */ 213*2f75a4aaSAlp Dener isf[n_isf]=low+i; ++n_isf; 214*2f75a4aaSAlp Dener } else if ((x[i] <= xl[i] + *bound_tol) && (g[i] > 0.0)) { 215*2f75a4aaSAlp Dener /* Lower bounded variables here */ 216*2f75a4aaSAlp Dener isl[n_isl]=low+i; ++n_isl; 217*2f75a4aaSAlp Dener } else if ((x[i] >= xu[i] - *bound_tol) && (g[i] < 0.0)) { 218*2f75a4aaSAlp Dener /* Upper bounded variables here */ 219*2f75a4aaSAlp Dener isu[n_isu]=low+i; ++n_isu; 220*2f75a4aaSAlp Dener } 221*2f75a4aaSAlp Dener } 222*2f75a4aaSAlp Dener 223*2f75a4aaSAlp Dener ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 224*2f75a4aaSAlp Dener ierr = VecRestoreArrayRead(XL, &xl);CHKERRQ(ierr); 225*2f75a4aaSAlp Dener ierr = VecRestoreArrayRead(XU, &xu);CHKERRQ(ierr); 226*2f75a4aaSAlp Dener ierr = VecRestoreArrayRead(G, &g);CHKERRQ(ierr); 227*2f75a4aaSAlp Dener } 228*2f75a4aaSAlp Dener 229*2f75a4aaSAlp Dener /* Create index set for lower bounded variables */ 230*2f75a4aaSAlp Dener ierr = ISDestroy(active_lower);CHKERRQ(ierr); 231*2f75a4aaSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isl, isl, PETSC_OWN_POINTER, active_lower);CHKERRQ(ierr); 232*2f75a4aaSAlp Dener /* Create index set for upper bounded variables */ 233*2f75a4aaSAlp Dener ierr = ISDestroy(active_upper);CHKERRQ(ierr); 234*2f75a4aaSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isu, isu, PETSC_OWN_POINTER, active_upper);CHKERRQ(ierr); 235*2f75a4aaSAlp Dener /* Create index set for fixed variables */ 236*2f75a4aaSAlp Dener ierr = ISDestroy(active_fixed);CHKERRQ(ierr); 237*2f75a4aaSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isf, isf, PETSC_OWN_POINTER, active_fixed);CHKERRQ(ierr); 238*2f75a4aaSAlp Dener 239*2f75a4aaSAlp Dener /* Create the combined active set */ 240*2f75a4aaSAlp Dener ierr = ISDestroy(active);CHKERRQ(ierr); 241*2f75a4aaSAlp Dener if (*active_lower && *active_upper && *active_fixed) { 242*2f75a4aaSAlp Dener /* All three types of active variables exist */ 243*2f75a4aaSAlp Dener const IS islist[3] = {*active_lower, *active_upper, *active_fixed}; 244*2f75a4aaSAlp Dener ierr = ISConcatenate(PetscObjectComm((PetscObject)X), 3, islist, active);CHKERRQ(ierr); 245*2f75a4aaSAlp Dener ierr = ISSort(*active);CHKERRQ(ierr); 246*2f75a4aaSAlp Dener } else if (*active_lower && *active_upper) { 247*2f75a4aaSAlp Dener /* Only lower and upper bounded active variables exist */ 248*2f75a4aaSAlp Dener ierr = ISSum(*active_lower, *active_upper, active);CHKERRQ(ierr); 249*2f75a4aaSAlp Dener } else if (*active_lower && *active_fixed) { 250*2f75a4aaSAlp Dener /* Only lower bounded and fixed active variables exist */ 251*2f75a4aaSAlp Dener ierr = ISSum(*active_lower, *active_fixed, active);CHKERRQ(ierr); 252*2f75a4aaSAlp Dener } else if (*active_upper && *active_fixed) { 253*2f75a4aaSAlp Dener /* Only upper bounded and fixed active variables exist */ 254*2f75a4aaSAlp Dener ierr = ISSum(*active_upper, *active_fixed, active);CHKERRQ(ierr); 255*2f75a4aaSAlp Dener } else if (*active_lower) { 256*2f75a4aaSAlp Dener /* Only lower bounded active variables exist */ 257*2f75a4aaSAlp Dener *active = *active_lower; 258*2f75a4aaSAlp Dener } else if (*active_upper) { 259*2f75a4aaSAlp Dener /* Only upper bounded active variables exist */ 260*2f75a4aaSAlp Dener *active = *active_upper; 261*2f75a4aaSAlp Dener } else if (*active_fixed) { 262*2f75a4aaSAlp Dener /* Only fixed active variables exist */ 263*2f75a4aaSAlp Dener *active = *active_fixed; 264*2f75a4aaSAlp Dener } 265*2f75a4aaSAlp Dener /* Create the inactive set */ 266*2f75a4aaSAlp Dener ierr = ISDestroy(inactive);CHKERRQ(ierr); 267*2f75a4aaSAlp Dener if (*active) { ierr = ISComplementVec(*active, X, inactive);CHKERRQ(ierr); } 268*2f75a4aaSAlp Dener 269*2f75a4aaSAlp Dener /* Clean up and exit */ 270*2f75a4aaSAlp Dener ierr = VecDestroy(&W);CHKERRQ(ierr); 271*2f75a4aaSAlp Dener PetscFunctionReturn(0); 272*2f75a4aaSAlp Dener } 273*2f75a4aaSAlp Dener 274*2f75a4aaSAlp Dener /*@C 275*2f75a4aaSAlp Dener TaoBoundStep - Ensures the correct zero or adjusted step direction 276*2f75a4aaSAlp Dener values for active variables. 277*2f75a4aaSAlp Dener 278*2f75a4aaSAlp Dener Input Parameters: 279*2f75a4aaSAlp Dener + X - solution vector 280*2f75a4aaSAlp Dener . XL - lower bound vector 281*2f75a4aaSAlp Dener . XU - upper bound vector 282*2f75a4aaSAlp Dener . active_lower - index set for lower bounded active variables 283*2f75a4aaSAlp Dener . active_upper - index set for lower bounded active variables 284*2f75a4aaSAlp Dener - active_fixed - index set for fixed active variables 285*2f75a4aaSAlp Dener 286*2f75a4aaSAlp Dener Output Parameters: 287*2f75a4aaSAlp Dener . S - step direction to be modified 288*2f75a4aaSAlp Dener @*/ 289*2f75a4aaSAlp Dener PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, Vec S) 290*2f75a4aaSAlp Dener { 291*2f75a4aaSAlp Dener PetscErrorCode ierr; 292*2f75a4aaSAlp Dener 293*2f75a4aaSAlp Dener Vec step_lower, step_upper, step_fixed; 294*2f75a4aaSAlp Dener Vec x_lower, x_upper; 295*2f75a4aaSAlp Dener Vec bound_lower, bound_upper; 296*2f75a4aaSAlp Dener 297*2f75a4aaSAlp Dener PetscFunctionBegin; 298*2f75a4aaSAlp Dener /* Adjust step for variables at the estimated lower bound */ 299*2f75a4aaSAlp Dener if (active_lower) { 300*2f75a4aaSAlp Dener ierr = VecGetSubVector(S, active_lower, &step_lower);CHKERRQ(ierr); 301*2f75a4aaSAlp Dener ierr = VecGetSubVector(X, active_lower, &x_lower);CHKERRQ(ierr); 302*2f75a4aaSAlp Dener ierr = VecGetSubVector(XL, active_lower, &bound_lower);CHKERRQ(ierr); 303*2f75a4aaSAlp Dener ierr = VecCopy(bound_lower, step_lower);CHKERRQ(ierr); 304*2f75a4aaSAlp Dener ierr = VecAXPY(step_lower, -1.0, x_lower);CHKERRQ(ierr); 305*2f75a4aaSAlp Dener ierr = VecRestoreSubVector(S, active_lower, &step_lower);CHKERRQ(ierr); 306*2f75a4aaSAlp Dener ierr = VecRestoreSubVector(X, active_lower, &x_lower);CHKERRQ(ierr); 307*2f75a4aaSAlp Dener ierr = VecRestoreSubVector(XL, active_lower, &bound_lower);CHKERRQ(ierr); 308*2f75a4aaSAlp Dener } 309*2f75a4aaSAlp Dener 310*2f75a4aaSAlp Dener /* Adjust step for the variables at the estimated upper bound */ 311*2f75a4aaSAlp Dener if (active_upper) { 312*2f75a4aaSAlp Dener ierr = VecGetSubVector(S, active_upper, &step_upper);CHKERRQ(ierr); 313*2f75a4aaSAlp Dener ierr = VecGetSubVector(X, active_upper, &x_upper);CHKERRQ(ierr); 314*2f75a4aaSAlp Dener ierr = VecGetSubVector(XU, active_upper, &bound_upper);CHKERRQ(ierr); 315*2f75a4aaSAlp Dener ierr = VecCopy(bound_upper, step_upper);CHKERRQ(ierr); 316*2f75a4aaSAlp Dener ierr = VecAXPY(step_upper, -1.0, x_upper);CHKERRQ(ierr); 317*2f75a4aaSAlp Dener ierr = VecRestoreSubVector(S, active_upper, &step_upper);CHKERRQ(ierr); 318*2f75a4aaSAlp Dener ierr = VecRestoreSubVector(X, active_upper, &x_upper);CHKERRQ(ierr); 319*2f75a4aaSAlp Dener ierr = VecRestoreSubVector(XU, active_upper, &bound_upper);CHKERRQ(ierr); 320*2f75a4aaSAlp Dener } 321*2f75a4aaSAlp Dener 322*2f75a4aaSAlp Dener /* Zero out step for fixed variables */ 323*2f75a4aaSAlp Dener if (active_fixed) { 324*2f75a4aaSAlp Dener ierr = VecGetSubVector(S, active_fixed, &step_fixed);CHKERRQ(ierr); 325*2f75a4aaSAlp Dener ierr = VecSet(step_fixed, 0.0);CHKERRQ(ierr); 326*2f75a4aaSAlp Dener ierr = VecRestoreSubVector(S, active_fixed, &step_fixed);CHKERRQ(ierr); 327*2f75a4aaSAlp Dener } 328*2f75a4aaSAlp Dener PetscFunctionReturn(0); 329*2f75a4aaSAlp Dener } 330