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 15a7e14dcfSSatish Balay Output Parameters: 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 21a7e14dcfSSatish Balay @*/ 223a831ad5SBarry Smith PetscErrorCode TaoVecGetSubVec(Vec vfull, IS is, TaoSubsetType reduced_type, PetscReal maskvalue, Vec *vreduced) 23a7e14dcfSSatish Balay { 24a7e14dcfSSatish Balay PetscErrorCode ierr; 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 38a7e14dcfSSatish Balay ierr = VecGetSize(vfull, &nfull);CHKERRQ(ierr); 39a7e14dcfSSatish Balay ierr = ISGetSize(is, &nreduced);CHKERRQ(ierr); 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay if (nreduced == nfull) { 42a7e14dcfSSatish Balay ierr = VecDestroy(vreduced);CHKERRQ(ierr); 43a7e14dcfSSatish Balay ierr = VecDuplicate(vfull,vreduced);CHKERRQ(ierr); 44a7e14dcfSSatish Balay ierr = VecCopy(vfull,*vreduced);CHKERRQ(ierr); 45a7e14dcfSSatish Balay } else { 46a7e14dcfSSatish Balay switch (reduced_type) { 47a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 48a7e14dcfSSatish Balay ierr = VecGetType(vfull,&vtype);CHKERRQ(ierr); 49a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(vfull,&flow,&fhigh);CHKERRQ(ierr); 50a7e14dcfSSatish Balay ierr = ISGetLocalSize(is,&nreduced_local);CHKERRQ(ierr); 51a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)vfull,&comm);CHKERRQ(ierr); 52a7e14dcfSSatish Balay if (*vreduced) { 53a7e14dcfSSatish Balay ierr = VecDestroy(vreduced);CHKERRQ(ierr); 54a7e14dcfSSatish Balay } 55a7e14dcfSSatish Balay ierr = VecCreate(comm,vreduced);CHKERRQ(ierr); 56a7e14dcfSSatish Balay ierr = VecSetType(*vreduced,vtype);CHKERRQ(ierr); 57a7e14dcfSSatish Balay 58a7e14dcfSSatish Balay ierr = VecSetSizes(*vreduced,nreduced_local,nreduced);CHKERRQ(ierr); 59a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(*vreduced,&rlow,&rhigh);CHKERRQ(ierr); 60a7e14dcfSSatish Balay ierr = ISCreateStride(comm,nreduced_local,rlow,1,&ident);CHKERRQ(ierr); 61a7e14dcfSSatish Balay ierr = VecScatterCreate(vfull,is,*vreduced,ident,&scatter);CHKERRQ(ierr); 62a7e14dcfSSatish Balay ierr = VecScatterBegin(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 63a7e14dcfSSatish Balay ierr = VecScatterEnd(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 64a7e14dcfSSatish Balay ierr = VecScatterDestroy(&scatter);CHKERRQ(ierr); 65a7e14dcfSSatish Balay ierr = ISDestroy(&ident);CHKERRQ(ierr); 66a7e14dcfSSatish Balay break; 67a7e14dcfSSatish Balay 68a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 69a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 70a7e14dcfSSatish Balay /* vr[i] = vf[i] if i in is 71a7e14dcfSSatish Balay vr[i] = 0 otherwise */ 726c4ed002SBarry Smith if (!*vreduced) { 73a7e14dcfSSatish Balay ierr = VecDuplicate(vfull,vreduced);CHKERRQ(ierr); 74a7e14dcfSSatish Balay } 75a7e14dcfSSatish Balay 76a7e14dcfSSatish Balay ierr = VecSet(*vreduced,maskvalue);CHKERRQ(ierr); 77a7e14dcfSSatish Balay ierr = ISGetLocalSize(is,&nlocal);CHKERRQ(ierr); 78a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(vfull,&flow,&fhigh);CHKERRQ(ierr); 79a7e14dcfSSatish Balay ierr = VecGetArray(vfull,&fv);CHKERRQ(ierr); 80a7e14dcfSSatish Balay ierr = VecGetArray(*vreduced,&rv);CHKERRQ(ierr); 81a7e14dcfSSatish Balay ierr = ISGetIndices(is,&s);CHKERRQ(ierr); 8253506e15SBarry Smith if (nlocal > (fhigh-flow)) SETERRQ2(PETSC_COMM_WORLD,1,"IS local size %d > Vec local size %d",nlocal,fhigh-flow); 83c4b75bccSAlp Dener for (i=0;i<nlocal;++i) { 84a7e14dcfSSatish Balay rv[s[i]-flow] = fv[s[i]-flow]; 85a7e14dcfSSatish Balay } 86a7e14dcfSSatish Balay ierr = ISRestoreIndices(is,&s);CHKERRQ(ierr); 87a7e14dcfSSatish Balay ierr = VecRestoreArray(vfull,&fv);CHKERRQ(ierr); 88a7e14dcfSSatish Balay ierr = VecRestoreArray(*vreduced,&rv);CHKERRQ(ierr); 89a7e14dcfSSatish Balay break; 90a7e14dcfSSatish Balay } 91a7e14dcfSSatish Balay } 92a7e14dcfSSatish Balay PetscFunctionReturn(0); 93a7e14dcfSSatish Balay } 94a7e14dcfSSatish Balay 95b98f30f2SJason Sarich /*@C 96b98f30f2SJason Sarich TaoMatGetSubMat - Gets a submatrix using the IS 97a7e14dcfSSatish Balay 98a7e14dcfSSatish Balay Input Parameters: 99a7e14dcfSSatish Balay + M - the full matrix (n x n) 100a7e14dcfSSatish Balay . is - the index set for the submatrix (both row and column index sets need to be the same) 101a7e14dcfSSatish Balay . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option 102a7e14dcfSSatish Balay - subset_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK, 103a7e14dcfSSatish Balay TAO_SUBSET_MATRIXFREE) 104a7e14dcfSSatish Balay 105a7e14dcfSSatish Balay Output Parameters: 106a7e14dcfSSatish Balay . Msub - the submatrix 107a7e14dcfSSatish Balay @*/ 108b98f30f2SJason Sarich PetscErrorCode TaoMatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub) 109a7e14dcfSSatish Balay { 110a7e14dcfSSatish Balay PetscErrorCode ierr; 111a7e14dcfSSatish Balay IS iscomp; 11262675beeSAlp Dener PetscBool flg = PETSC_TRUE; 11353506e15SBarry Smith 114a7e14dcfSSatish Balay PetscFunctionBegin; 115a7e14dcfSSatish Balay PetscValidHeaderSpecific(M,MAT_CLASSID,1); 116a7e14dcfSSatish Balay PetscValidHeaderSpecific(is,IS_CLASSID,2); 117a7e14dcfSSatish Balay ierr = MatDestroy(Msub);CHKERRQ(ierr); 118a7e14dcfSSatish Balay switch (subset_type) { 119a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 1207dae84e0SHong Zhang ierr = MatCreateSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub);CHKERRQ(ierr); 121a7e14dcfSSatish Balay break; 122a7e14dcfSSatish Balay 123a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 124a7e14dcfSSatish Balay /* Get Reduced Hessian 125a7e14dcfSSatish Balay Msub[i,j] = M[i,j] if i,j in Free_Local or i==j 126a7e14dcfSSatish Balay Msub[i,j] = 0 if i!=j and i or j not in Free_Local 127a7e14dcfSSatish Balay */ 1281a1499c8SBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)M),NULL,NULL,NULL);CHKERRQ(ierr); 12962675beeSAlp Dener ierr = PetscOptionsBool("-overwrite_hessian","modify the existing hessian matrix when computing submatrices","TaoSubsetType",flg,&flg,NULL);CHKERRQ(ierr); 130302440fdSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 1318afaa268SBarry Smith if (flg) { 132a7e14dcfSSatish Balay ierr = MatDuplicate(M, MAT_COPY_VALUES, Msub);CHKERRQ(ierr); 133a7e14dcfSSatish Balay } else { 134a7e14dcfSSatish Balay /* Act on hessian directly (default) */ 135a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)M);CHKERRQ(ierr); 136a7e14dcfSSatish Balay *Msub = M; 137a7e14dcfSSatish Balay } 138a7e14dcfSSatish Balay /* Save the diagonal to temporary vector */ 139a7e14dcfSSatish Balay ierr = MatGetDiagonal(*Msub,v1);CHKERRQ(ierr); 140a7e14dcfSSatish Balay 141a7e14dcfSSatish Balay /* Zero out rows and columns */ 1424473680cSBarry Smith ierr = ISComplementVec(is,v1,&iscomp);CHKERRQ(ierr); 143a7e14dcfSSatish Balay 144a7e14dcfSSatish Balay /* Use v1 instead of 0 here because of PETSc bug */ 145a7e14dcfSSatish Balay ierr = MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1);CHKERRQ(ierr); 146a7e14dcfSSatish Balay 147a7e14dcfSSatish Balay ierr = ISDestroy(&iscomp);CHKERRQ(ierr); 148a7e14dcfSSatish Balay break; 149a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 1504473680cSBarry Smith ierr = ISComplementVec(is,v1,&iscomp);CHKERRQ(ierr); 151a7e14dcfSSatish Balay ierr = MatCreateSubMatrixFree(M,iscomp,iscomp,Msub);CHKERRQ(ierr); 152a7e14dcfSSatish Balay ierr = ISDestroy(&iscomp);CHKERRQ(ierr); 153a7e14dcfSSatish Balay break; 154a7e14dcfSSatish Balay } 155a7e14dcfSSatish Balay PetscFunctionReturn(0); 156a7e14dcfSSatish Balay } 1572f75a4aaSAlp Dener 1582f75a4aaSAlp Dener /*@C 1592f75a4aaSAlp Dener TaoEstimateActiveBounds - Generates index sets for variables at the lower and upper 1602f75a4aaSAlp Dener bounds, as well as fixed variables where lower and upper bounds equal each other. 1612f75a4aaSAlp Dener 1622f75a4aaSAlp Dener Input Parameters: 1632f75a4aaSAlp Dener + X - solution vector 1642f75a4aaSAlp Dener . XL - lower bound vector 1652f75a4aaSAlp Dener . XU - upper bound vector 1662f75a4aaSAlp Dener . G - unprojected gradient 1670a4511e9SAlp Dener . S - step direction with which the active bounds will be estimated 1680a4511e9SAlp Dener - steplen - the step length at which the active bounds will be estimated (needs to be conservative) 1692f75a4aaSAlp Dener 1702f75a4aaSAlp Dener Output Parameters: 1712f75a4aaSAlp Dener . bound_tol - tolerance for for the bound estimation 1722f75a4aaSAlp Dener . active_lower - index set for active variables at the lower bound 1732f75a4aaSAlp Dener . active_upper - index set for active variables at the upper bound 1742f75a4aaSAlp Dener . active_fixed - index set for fixed variables 1752f75a4aaSAlp Dener . active - index set for all active variables 1762f75a4aaSAlp Dener . inactive - complementary index set for inactive variables 1770a4511e9SAlp Dener 1780a4511e9SAlp Dener Notes: 1790a4511e9SAlp Dener This estimation is based on Bertsekas' method, with a built in diagonal scaling value of 1.0e-3. 1800a4511e9SAlp Dener 1812f75a4aaSAlp Dener @*/ 182c4b75bccSAlp Dener PetscErrorCode TaoEstimateActiveBounds(Vec X, Vec XL, Vec XU, Vec G, Vec S, Vec W, PetscReal steplen, PetscReal *bound_tol, 1832f75a4aaSAlp Dener IS *active_lower, IS *active_upper, IS *active_fixed, IS *active, IS *inactive) 1842f75a4aaSAlp Dener { 1852f75a4aaSAlp Dener PetscErrorCode ierr; 1862f75a4aaSAlp Dener 1872f75a4aaSAlp Dener PetscReal wnorm; 188c4b75bccSAlp Dener PetscReal mach_eps = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); 189c4b75bccSAlp Dener PetscInt i, n_isl=0, n_isu=0, n_isf=0, n_isa=0, n_isi=0; 190c4b75bccSAlp Dener PetscInt N_isl, N_isu, N_isf, N_isa, N_isi; 191c4b75bccSAlp Dener PetscInt n, low, high, nDiff; 192c4b75bccSAlp Dener PetscInt *isl=NULL, *isu=NULL, *isf=NULL, *isa=NULL, *isi=NULL; 1932f75a4aaSAlp Dener const PetscScalar *xl, *xu, *x, *g; 1942f75a4aaSAlp Dener 1952f75a4aaSAlp Dener PetscFunctionBegin; 196c4b75bccSAlp Dener PetscValidHeaderSpecific(X,VEC_CLASSID,1); 197c4b75bccSAlp Dener PetscValidHeaderSpecific(XL,VEC_CLASSID,2); 198c4b75bccSAlp Dener PetscValidHeaderSpecific(XU,VEC_CLASSID,3); 199c4b75bccSAlp Dener PetscValidHeaderSpecific(G,VEC_CLASSID,4); 200c4b75bccSAlp Dener PetscValidHeaderSpecific(S,VEC_CLASSID,5); 201c4b75bccSAlp Dener PetscValidHeaderSpecific(W,VEC_CLASSID,6); 202c4b75bccSAlp Dener 203c4b75bccSAlp Dener PetscValidType(X,1); 204c4b75bccSAlp Dener PetscValidType(XL,2); 205c4b75bccSAlp Dener PetscValidType(XU,3); 206c4b75bccSAlp Dener PetscValidType(G,4); 207c4b75bccSAlp Dener PetscValidType(S,5); 208c4b75bccSAlp Dener PetscValidType(W,6); 209c4b75bccSAlp Dener PetscCheckSameType(X,1,XL,2); 210c4b75bccSAlp Dener PetscCheckSameType(X,1,XU,3); 211c4b75bccSAlp Dener PetscCheckSameType(X,1,G,4); 212c4b75bccSAlp Dener PetscCheckSameType(X,1,S,5); 213c4b75bccSAlp Dener PetscCheckSameType(X,1,W,6); 214c4b75bccSAlp Dener PetscCheckSameComm(X,1,XL,2); 215c4b75bccSAlp Dener PetscCheckSameComm(X,1,XU,3); 216c4b75bccSAlp Dener PetscCheckSameComm(X,1,G,4); 217c4b75bccSAlp Dener PetscCheckSameComm(X,1,S,5); 218c4b75bccSAlp Dener PetscCheckSameComm(X,1,W,6); 219c4b75bccSAlp Dener VecCheckSameSize(X,1,XL,2); 220c4b75bccSAlp Dener VecCheckSameSize(X,1,XU,3); 221c4b75bccSAlp Dener VecCheckSameSize(X,1,G,4); 222c4b75bccSAlp Dener VecCheckSameSize(X,1,S,5); 223c4b75bccSAlp Dener VecCheckSameSize(X,1,W,6); 224c4b75bccSAlp Dener 2252f75a4aaSAlp Dener /* Update the tolerance for bound detection (this is based on Bertsekas' method) */ 226c4b75bccSAlp Dener ierr = VecCopy(X, W);CHKERRQ(ierr); 227c4b75bccSAlp Dener ierr = VecAXPBY(W, steplen, 1.0, S);CHKERRQ(ierr); 228c4b75bccSAlp Dener ierr = TaoBoundSolution(XL, XU, W, &nDiff);CHKERRQ(ierr); 2292f75a4aaSAlp Dener ierr = VecAXPBY(W, 1.0, -1.0, X);CHKERRQ(ierr); 2302f75a4aaSAlp Dener ierr = VecNorm(W, NORM_2, &wnorm);CHKERRQ(ierr); 2312f75a4aaSAlp Dener *bound_tol = PetscMin(*bound_tol, wnorm); 2322f75a4aaSAlp Dener 2332f75a4aaSAlp Dener ierr = VecGetOwnershipRange(X, &low, &high);CHKERRQ(ierr); 2342f75a4aaSAlp Dener ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr); 2352f75a4aaSAlp Dener if (n>0){ 2362f75a4aaSAlp Dener ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 2372f75a4aaSAlp Dener ierr = VecGetArrayRead(XL, &xl);CHKERRQ(ierr); 2382f75a4aaSAlp Dener ierr = VecGetArrayRead(XU, &xu);CHKERRQ(ierr); 2392f75a4aaSAlp Dener ierr = VecGetArrayRead(G, &g);CHKERRQ(ierr); 2402f75a4aaSAlp Dener 2412f75a4aaSAlp Dener /* Loop over variables and categorize the indexes */ 2422f75a4aaSAlp Dener ierr = PetscMalloc1(n, &isl);CHKERRQ(ierr); 2432f75a4aaSAlp Dener ierr = PetscMalloc1(n, &isu);CHKERRQ(ierr); 2442f75a4aaSAlp Dener ierr = PetscMalloc1(n, &isf);CHKERRQ(ierr); 245c4b75bccSAlp Dener ierr = PetscMalloc1(n, &isa);CHKERRQ(ierr); 246c4b75bccSAlp Dener ierr = PetscMalloc1(n, &isi);CHKERRQ(ierr); 247c4b75bccSAlp Dener for (i=0; i<n; ++i) { 248*66b4d56fSTodd Munson if (xl[i] == xu[i]) { 249c4b75bccSAlp Dener /* Fixed variables */ 2502f75a4aaSAlp Dener isf[n_isf]=low+i; ++n_isf; 251c4b75bccSAlp Dener isa[n_isa]=low+i; ++n_isa; 252*66b4d56fSTodd Munson } else if ((xl[i] > PETSC_NINFINITY) && (x[i] <= xl[i] + *bound_tol) && (g[i] >= -mach_eps)) { 253c4b75bccSAlp Dener /* Lower bounded variables */ 2542f75a4aaSAlp Dener isl[n_isl]=low+i; ++n_isl; 255c4b75bccSAlp Dener isa[n_isa]=low+i; ++n_isa; 256*66b4d56fSTodd Munson } else if ((xu[i] < PETSC_INFINITY) && (x[i] >= xu[i] - *bound_tol) && (g[i] <= mach_eps)) { 257c4b75bccSAlp Dener /* Upper bounded variables */ 2582f75a4aaSAlp Dener isu[n_isu]=low+i; ++n_isu; 259c4b75bccSAlp Dener isa[n_isa]=low+i; ++n_isa; 260c4b75bccSAlp Dener } else { 261c4b75bccSAlp Dener /* Inactive variables */ 262c4b75bccSAlp Dener isi[n_isi]=low+i; ++n_isi; 2632f75a4aaSAlp Dener } 2642f75a4aaSAlp Dener } 2652f75a4aaSAlp Dener 2662f75a4aaSAlp Dener ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 2672f75a4aaSAlp Dener ierr = VecRestoreArrayRead(XL, &xl);CHKERRQ(ierr); 2682f75a4aaSAlp Dener ierr = VecRestoreArrayRead(XU, &xu);CHKERRQ(ierr); 2692f75a4aaSAlp Dener ierr = VecRestoreArrayRead(G, &g);CHKERRQ(ierr); 2702f75a4aaSAlp Dener } 2712f75a4aaSAlp Dener 272c4b75bccSAlp Dener /* Clear all index sets */ 2732f75a4aaSAlp Dener ierr = ISDestroy(active_lower);CHKERRQ(ierr); 2742f75a4aaSAlp Dener ierr = ISDestroy(active_upper);CHKERRQ(ierr); 2752f75a4aaSAlp Dener ierr = ISDestroy(active_fixed);CHKERRQ(ierr); 2762f75a4aaSAlp Dener ierr = ISDestroy(active);CHKERRQ(ierr); 2772f75a4aaSAlp Dener ierr = ISDestroy(inactive);CHKERRQ(ierr); 278c4b75bccSAlp Dener 279c4b75bccSAlp Dener /* Collect global sizes */ 280c4b75bccSAlp Dener ierr = MPIU_Allreduce(&n_isl, &N_isl, 1, MPI_INT, MPI_SUM, PetscObjectComm((PetscObject)X));CHKERRQ(ierr); 281c4b75bccSAlp Dener ierr = MPIU_Allreduce(&n_isu, &N_isu, 1, MPI_INT, MPI_SUM, PetscObjectComm((PetscObject)X));CHKERRQ(ierr); 282c4b75bccSAlp Dener ierr = MPIU_Allreduce(&n_isf, &N_isf, 1, MPI_INT, MPI_SUM, PetscObjectComm((PetscObject)X));CHKERRQ(ierr); 283c4b75bccSAlp Dener ierr = MPIU_Allreduce(&n_isa, &N_isa, 1, MPI_INT, MPI_SUM, PetscObjectComm((PetscObject)X));CHKERRQ(ierr); 284c4b75bccSAlp Dener ierr = MPIU_Allreduce(&n_isi, &N_isi, 1, MPI_INT, MPI_SUM, PetscObjectComm((PetscObject)X));CHKERRQ(ierr); 285c4b75bccSAlp Dener 286c4b75bccSAlp Dener /* Create index set for lower bounded variables */ 287c4b75bccSAlp Dener if (N_isl > 0) { 288c4b75bccSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isl, isl, PETSC_OWN_POINTER, active_lower);CHKERRQ(ierr); 289c4b75bccSAlp Dener } 290c4b75bccSAlp Dener /* Create index set for upper bounded variables */ 291c4b75bccSAlp Dener if (N_isu > 0) { 292c4b75bccSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isu, isu, PETSC_OWN_POINTER, active_upper);CHKERRQ(ierr); 293c4b75bccSAlp Dener } 294c4b75bccSAlp Dener /* Create index set for fixed variables */ 295c4b75bccSAlp Dener if (N_isf > 0) { 296c4b75bccSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isf, isf, PETSC_OWN_POINTER, active_fixed);CHKERRQ(ierr); 297c4b75bccSAlp Dener } 298c4b75bccSAlp Dener /* Create index set for all actively bounded variables */ 299c4b75bccSAlp Dener if (N_isa > 0) { 300c4b75bccSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isa, isa, PETSC_OWN_POINTER, active);CHKERRQ(ierr); 301c4b75bccSAlp Dener } 302c4b75bccSAlp Dener /* Create index set for all inactive variables */ 303c4b75bccSAlp Dener if (N_isi > 0) { 304c4b75bccSAlp Dener ierr = ISCreateGeneral(PetscObjectComm((PetscObject)X), n_isi, isi, PETSC_OWN_POINTER, inactive);CHKERRQ(ierr); 305c4b75bccSAlp Dener } 3062f75a4aaSAlp Dener 3072f75a4aaSAlp Dener /* Clean up and exit */ 3082f75a4aaSAlp Dener PetscFunctionReturn(0); 3092f75a4aaSAlp Dener } 3102f75a4aaSAlp Dener 3112f75a4aaSAlp Dener /*@C 3122f75a4aaSAlp Dener TaoBoundStep - Ensures the correct zero or adjusted step direction 3132f75a4aaSAlp Dener values for active variables. 3142f75a4aaSAlp Dener 3152f75a4aaSAlp Dener Input Parameters: 3162f75a4aaSAlp Dener + X - solution vector 3172f75a4aaSAlp Dener . XL - lower bound vector 3182f75a4aaSAlp Dener . XU - upper bound vector 3192f75a4aaSAlp Dener . active_lower - index set for lower bounded active variables 3202f75a4aaSAlp Dener . active_upper - index set for lower bounded active variables 3212f75a4aaSAlp Dener - active_fixed - index set for fixed active variables 3222f75a4aaSAlp Dener 3232f75a4aaSAlp Dener Output Parameters: 3242f75a4aaSAlp Dener . S - step direction to be modified 3252f75a4aaSAlp Dener @*/ 326c4b75bccSAlp Dener PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, PetscReal scale, Vec S) 3272f75a4aaSAlp Dener { 3282f75a4aaSAlp Dener PetscErrorCode ierr; 3292f75a4aaSAlp Dener 3302f75a4aaSAlp Dener Vec step_lower, step_upper, step_fixed; 3312f75a4aaSAlp Dener Vec x_lower, x_upper; 3322f75a4aaSAlp Dener Vec bound_lower, bound_upper; 3332f75a4aaSAlp Dener 3342f75a4aaSAlp Dener PetscFunctionBegin; 3352f75a4aaSAlp Dener /* Adjust step for variables at the estimated lower bound */ 3362f75a4aaSAlp Dener if (active_lower) { 3372f75a4aaSAlp Dener ierr = VecGetSubVector(S, active_lower, &step_lower);CHKERRQ(ierr); 3382f75a4aaSAlp Dener ierr = VecGetSubVector(X, active_lower, &x_lower);CHKERRQ(ierr); 3392f75a4aaSAlp Dener ierr = VecGetSubVector(XL, active_lower, &bound_lower);CHKERRQ(ierr); 3402f75a4aaSAlp Dener ierr = VecCopy(bound_lower, step_lower);CHKERRQ(ierr); 3412f75a4aaSAlp Dener ierr = VecAXPY(step_lower, -1.0, x_lower);CHKERRQ(ierr); 342c4b75bccSAlp Dener ierr = VecScale(step_lower, scale);CHKERRQ(ierr); 3432f75a4aaSAlp Dener ierr = VecRestoreSubVector(S, active_lower, &step_lower);CHKERRQ(ierr); 3442f75a4aaSAlp Dener ierr = VecRestoreSubVector(X, active_lower, &x_lower);CHKERRQ(ierr); 3452f75a4aaSAlp Dener ierr = VecRestoreSubVector(XL, active_lower, &bound_lower);CHKERRQ(ierr); 3462f75a4aaSAlp Dener } 3472f75a4aaSAlp Dener 3482f75a4aaSAlp Dener /* Adjust step for the variables at the estimated upper bound */ 3492f75a4aaSAlp Dener if (active_upper) { 3502f75a4aaSAlp Dener ierr = VecGetSubVector(S, active_upper, &step_upper);CHKERRQ(ierr); 3512f75a4aaSAlp Dener ierr = VecGetSubVector(X, active_upper, &x_upper);CHKERRQ(ierr); 3522f75a4aaSAlp Dener ierr = VecGetSubVector(XU, active_upper, &bound_upper);CHKERRQ(ierr); 3532f75a4aaSAlp Dener ierr = VecCopy(bound_upper, step_upper);CHKERRQ(ierr); 3542f75a4aaSAlp Dener ierr = VecAXPY(step_upper, -1.0, x_upper);CHKERRQ(ierr); 355c4b75bccSAlp Dener ierr = VecScale(step_upper, scale);CHKERRQ(ierr); 3562f75a4aaSAlp Dener ierr = VecRestoreSubVector(S, active_upper, &step_upper);CHKERRQ(ierr); 3572f75a4aaSAlp Dener ierr = VecRestoreSubVector(X, active_upper, &x_upper);CHKERRQ(ierr); 3582f75a4aaSAlp Dener ierr = VecRestoreSubVector(XU, active_upper, &bound_upper);CHKERRQ(ierr); 3592f75a4aaSAlp Dener } 3602f75a4aaSAlp Dener 3612f75a4aaSAlp Dener /* Zero out step for fixed variables */ 3622f75a4aaSAlp Dener if (active_fixed) { 3632f75a4aaSAlp Dener ierr = VecGetSubVector(S, active_fixed, &step_fixed);CHKERRQ(ierr); 3642f75a4aaSAlp Dener ierr = VecSet(step_fixed, 0.0);CHKERRQ(ierr); 3652f75a4aaSAlp Dener ierr = VecRestoreSubVector(S, active_fixed, &step_fixed);CHKERRQ(ierr); 3662f75a4aaSAlp Dener } 3672f75a4aaSAlp Dener PetscFunctionReturn(0); 3682f75a4aaSAlp Dener } 369c4b75bccSAlp Dener 370c4b75bccSAlp Dener /*@C 371c4b75bccSAlp Dener TaoBoundSolution - Ensures that the solution vector is snapped into the bounds. 372c4b75bccSAlp Dener 373c4b75bccSAlp Dener Input Parameters: 374c4b75bccSAlp Dener + XL - lower bound vector 375c4b75bccSAlp Dener . XU - upper bound vector 376c4b75bccSAlp Dener - X - solution vector 377c4b75bccSAlp Dener 378c4b75bccSAlp Dener Output Parameters: 379c4b75bccSAlp Dener . X - modified solution vector 380c4b75bccSAlp Dener @*/ 381c4b75bccSAlp Dener PetscErrorCode TaoBoundSolution(Vec XL, Vec XU, Vec X, PetscInt *nDiff) 382c4b75bccSAlp Dener { 383c4b75bccSAlp Dener PetscErrorCode ierr; 384c4b75bccSAlp Dener PetscInt i,n,low,high,nDiff_loc=0; 385c4b75bccSAlp Dener PetscScalar *x; 386c4b75bccSAlp Dener const PetscScalar *xl,*xu; 387c4b75bccSAlp Dener 388c4b75bccSAlp Dener PetscFunctionBegin; 389c4b75bccSAlp Dener PetscValidHeaderSpecific(X,VEC_CLASSID,1); 390c4b75bccSAlp Dener PetscValidHeaderSpecific(XL,VEC_CLASSID,2); 391c4b75bccSAlp Dener PetscValidHeaderSpecific(XU,VEC_CLASSID,3); 392c4b75bccSAlp Dener 393c4b75bccSAlp Dener PetscValidType(X,1); 394c4b75bccSAlp Dener PetscValidType(XL,2); 395c4b75bccSAlp Dener PetscValidType(XU,3); 396c4b75bccSAlp Dener PetscCheckSameType(X,1,XL,2); 397c4b75bccSAlp Dener PetscCheckSameType(X,1,XU,3); 398c4b75bccSAlp Dener PetscCheckSameComm(X,1,XL,2); 399c4b75bccSAlp Dener PetscCheckSameComm(X,1,XU,3); 400c4b75bccSAlp Dener VecCheckSameSize(X,1,XL,2); 401c4b75bccSAlp Dener VecCheckSameSize(X,1,XU,3); 402c4b75bccSAlp Dener 403c4b75bccSAlp Dener ierr = VecGetOwnershipRange(X,&low,&high);CHKERRQ(ierr); 404c4b75bccSAlp Dener ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 405c4b75bccSAlp Dener if (n>0){ 406c4b75bccSAlp Dener ierr = VecGetArray(X, &x); 407c4b75bccSAlp Dener ierr = VecGetArrayRead(XL, &xl); 408c4b75bccSAlp Dener ierr = VecGetArrayRead(XU, &xu); 409c4b75bccSAlp Dener 410c4b75bccSAlp Dener for (i=0;i<n;++i){ 411c4b75bccSAlp Dener if ((x[i] < xl[i]) && (xl[i] > PETSC_NINFINITY)) { 412c4b75bccSAlp Dener x[i] = xl[i]; ++nDiff_loc; 413c4b75bccSAlp Dener } else if ((x[i] > xu[i]) && (xu[i] < PETSC_INFINITY)) { 414c4b75bccSAlp Dener x[i] = xu[i]; ++nDiff_loc; 415c4b75bccSAlp Dener } 416c4b75bccSAlp Dener } 417c4b75bccSAlp Dener 418c4b75bccSAlp Dener ierr = VecRestoreArray(X, &x); 419c4b75bccSAlp Dener ierr = VecRestoreArrayRead(XL, &xl); 420c4b75bccSAlp Dener ierr = VecRestoreArrayRead(XU, &xu); 421c4b75bccSAlp Dener } 422c4b75bccSAlp Dener ierr = MPIU_Allreduce(&nDiff_loc, nDiff, 1, MPI_INT, MPI_SUM, PetscObjectComm((PetscObject)X));CHKERRQ(ierr); 423c4b75bccSAlp Dener PetscFunctionReturn(0); 424c4b75bccSAlp Dener } 425