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 @*/ 233a831ad5SBarry Smith PetscErrorCode TaoVecGetSubVec(Vec vfull, IS is, TaoSubsetType reduced_type, PetscReal maskvalue, Vec *vreduced) 24a7e14dcfSSatish Balay { 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)); 52a7e14dcfSSatish Balay if (*vreduced) { 539566063dSJacob Faibussowitsch PetscCall(VecDestroy(vreduced)); 54a7e14dcfSSatish Balay } 559566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,vreduced)); 569566063dSJacob Faibussowitsch PetscCall(VecSetType(*vreduced,vtype)); 57a7e14dcfSSatish Balay 589566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*vreduced,nreduced_local,nreduced)); 599566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(*vreduced,&rlow,&rhigh)); 609566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm,nreduced_local,rlow,1,&ident)); 619566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vfull,is,*vreduced,ident,&scatter)); 629566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD)); 639566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD)); 649566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scatter)); 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ident)); 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) { 739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vfull,vreduced)); 74a7e14dcfSSatish Balay } 75a7e14dcfSSatish Balay 769566063dSJacob Faibussowitsch PetscCall(VecSet(*vreduced,maskvalue)); 779566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is,&nlocal)); 789566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vfull,&flow,&fhigh)); 799566063dSJacob Faibussowitsch PetscCall(VecGetArray(vfull,&fv)); 809566063dSJacob Faibussowitsch PetscCall(VecGetArray(*vreduced,&rv)); 819566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is,&s)); 8263a3b9bcSJacob 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); 83c4b75bccSAlp Dener for (i=0;i<nlocal;++i) { 84a7e14dcfSSatish Balay rv[s[i]-flow] = fv[s[i]-flow]; 85a7e14dcfSSatish Balay } 869566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is,&s)); 879566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vfull,&fv)); 889566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*vreduced,&rv)); 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 102b6a6cedcSAlp Dener - subset_type <TAO_SUBSET_SUBVEC,TAO_SUBSET_MASK,TAO_SUBSET_MATRIXFREE> - the method TAO is using for subsetting 103a7e14dcfSSatish Balay 10497bb3fdcSJose E. Roman Output Parameter: 105a7e14dcfSSatish Balay . Msub - the submatrix 106b6a6cedcSAlp Dener 107b6a6cedcSAlp Dener Level: developer 108a7e14dcfSSatish Balay @*/ 109b98f30f2SJason Sarich PetscErrorCode TaoMatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub) 110a7e14dcfSSatish Balay { 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); 1179566063dSJacob Faibussowitsch PetscCall(MatDestroy(Msub)); 118a7e14dcfSSatish Balay switch (subset_type) { 119a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 1209566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub)); 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 */ 128d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)M); 1299566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-overwrite_hessian","modify the existing hessian matrix when computing submatrices","TaoSubsetType",flg,&flg,NULL)); 130d0609cedSBarry Smith PetscOptionsEnd(); 1318afaa268SBarry Smith if (flg) { 1329566063dSJacob Faibussowitsch PetscCall(MatDuplicate(M, MAT_COPY_VALUES, Msub)); 133a7e14dcfSSatish Balay } else { 134a7e14dcfSSatish Balay /* Act on hessian directly (default) */ 1359566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)M)); 136a7e14dcfSSatish Balay *Msub = M; 137a7e14dcfSSatish Balay } 138a7e14dcfSSatish Balay /* Save the diagonal to temporary vector */ 1399566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(*Msub,v1)); 140a7e14dcfSSatish Balay 141a7e14dcfSSatish Balay /* Zero out rows and columns */ 1429566063dSJacob Faibussowitsch PetscCall(ISComplementVec(is,v1,&iscomp)); 143a7e14dcfSSatish Balay 144a7e14dcfSSatish Balay /* Use v1 instead of 0 here because of PETSc bug */ 1459566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1)); 146a7e14dcfSSatish Balay 1479566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscomp)); 148a7e14dcfSSatish Balay break; 149a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 1509566063dSJacob Faibussowitsch PetscCall(ISComplementVec(is,v1,&iscomp)); 1519566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrixFree(M,iscomp,iscomp,Msub)); 1529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscomp)); 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 168b6a6cedcSAlp Dener . W - work vector of type and size of X 1690a4511e9SAlp Dener - steplen - the step length at which the active bounds will be estimated (needs to be conservative) 1702f75a4aaSAlp Dener 1712f75a4aaSAlp Dener Output Parameters: 172*76be6f4fSStefano Zampini + bound_tol - tolerance for the bound estimation 1732f75a4aaSAlp Dener . active_lower - index set for active variables at the lower bound 1742f75a4aaSAlp Dener . active_upper - index set for active variables at the upper bound 1752f75a4aaSAlp Dener . active_fixed - index set for fixed variables 1762f75a4aaSAlp Dener . active - index set for all active variables 177b6a6cedcSAlp Dener - inactive - complementary index set for inactive variables 1780a4511e9SAlp Dener 1790a4511e9SAlp Dener Notes: 1800a4511e9SAlp Dener This estimation is based on Bertsekas' method, with a built in diagonal scaling value of 1.0e-3. 1810a4511e9SAlp Dener 182b6a6cedcSAlp Dener Level: developer 1832f75a4aaSAlp Dener @*/ 184c4b75bccSAlp Dener PetscErrorCode TaoEstimateActiveBounds(Vec X, Vec XL, Vec XU, Vec G, Vec S, Vec W, PetscReal steplen, PetscReal *bound_tol, 1852f75a4aaSAlp Dener IS *active_lower, IS *active_upper, IS *active_fixed, IS *active, IS *inactive) 1862f75a4aaSAlp Dener { 1872f75a4aaSAlp Dener PetscReal wnorm; 18889da521bSAlp Dener PetscReal zero = 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; 1946519dc0cSAlp Dener MPI_Comm comm = PetscObjectComm((PetscObject)X); 1952f75a4aaSAlp Dener 1962f75a4aaSAlp Dener PetscFunctionBegin; 197c4b75bccSAlp Dener PetscValidHeaderSpecific(X,VEC_CLASSID,1); 198*76be6f4fSStefano Zampini if (XL) PetscValidHeaderSpecific(XL,VEC_CLASSID,2); 199*76be6f4fSStefano Zampini if (XU) PetscValidHeaderSpecific(XU,VEC_CLASSID,3); 200c4b75bccSAlp Dener PetscValidHeaderSpecific(G,VEC_CLASSID,4); 201c4b75bccSAlp Dener PetscValidHeaderSpecific(S,VEC_CLASSID,5); 202c4b75bccSAlp Dener PetscValidHeaderSpecific(W,VEC_CLASSID,6); 203c4b75bccSAlp Dener 204*76be6f4fSStefano Zampini if (XL) PetscCheckSameType(X,1,XL,2); 205*76be6f4fSStefano Zampini if (XU) PetscCheckSameType(X,1,XU,3); 206c4b75bccSAlp Dener PetscCheckSameType(X,1,G,4); 207c4b75bccSAlp Dener PetscCheckSameType(X,1,S,5); 208c4b75bccSAlp Dener PetscCheckSameType(X,1,W,6); 209*76be6f4fSStefano Zampini if (XL) PetscCheckSameComm(X,1,XL,2); 210*76be6f4fSStefano Zampini if (XU) PetscCheckSameComm(X,1,XU,3); 211c4b75bccSAlp Dener PetscCheckSameComm(X,1,G,4); 212c4b75bccSAlp Dener PetscCheckSameComm(X,1,S,5); 213c4b75bccSAlp Dener PetscCheckSameComm(X,1,W,6); 214*76be6f4fSStefano Zampini if (XL) VecCheckSameSize(X,1,XL,2); 215*76be6f4fSStefano Zampini if (XU) VecCheckSameSize(X,1,XU,3); 216c4b75bccSAlp Dener VecCheckSameSize(X,1,G,4); 217c4b75bccSAlp Dener VecCheckSameSize(X,1,S,5); 218c4b75bccSAlp Dener VecCheckSameSize(X,1,W,6); 219c4b75bccSAlp Dener 2202f75a4aaSAlp Dener /* Update the tolerance for bound detection (this is based on Bertsekas' method) */ 2219566063dSJacob Faibussowitsch PetscCall(VecCopy(X, W)); 2229566063dSJacob Faibussowitsch PetscCall(VecAXPBY(W, steplen, 1.0, S)); 2239566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(W, XL, XU, 0.0, &nDiff, W)); 2249566063dSJacob Faibussowitsch PetscCall(VecAXPBY(W, 1.0, -1.0, X)); 2259566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 2262f75a4aaSAlp Dener *bound_tol = PetscMin(*bound_tol, wnorm); 2272f75a4aaSAlp Dener 228*76be6f4fSStefano Zampini /* Clear all index sets */ 229*76be6f4fSStefano Zampini PetscCall(ISDestroy(active_lower)); 230*76be6f4fSStefano Zampini PetscCall(ISDestroy(active_upper)); 231*76be6f4fSStefano Zampini PetscCall(ISDestroy(active_fixed)); 232*76be6f4fSStefano Zampini PetscCall(ISDestroy(active)); 233*76be6f4fSStefano Zampini PetscCall(ISDestroy(inactive)); 234*76be6f4fSStefano Zampini 2359566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &low, &high)); 2369566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 237*76be6f4fSStefano Zampini if (!XL && !XU) { 238*76be6f4fSStefano Zampini PetscCall(ISCreateStride(comm,n,low,1,inactive)); 239*76be6f4fSStefano Zampini PetscFunctionReturn(0); 240*76be6f4fSStefano Zampini } 2412f75a4aaSAlp Dener if (n>0) { 2429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2439566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XL, &xl)); 2449566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XU, &xu)); 2459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(G, &g)); 2462f75a4aaSAlp Dener 2472f75a4aaSAlp Dener /* Loop over variables and categorize the indexes */ 2489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isl)); 2499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isu)); 2509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isf)); 2519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isa)); 2529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &isi)); 253c4b75bccSAlp Dener for (i=0; i<n; ++i) { 25466b4d56fSTodd Munson if (xl[i] == xu[i]) { 255c4b75bccSAlp Dener /* Fixed variables */ 2562f75a4aaSAlp Dener isf[n_isf]=low+i; ++n_isf; 257c4b75bccSAlp Dener isa[n_isa]=low+i; ++n_isa; 258*76be6f4fSStefano Zampini } else if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + *bound_tol && g[i] > zero) { 259c4b75bccSAlp Dener /* Lower bounded variables */ 2602f75a4aaSAlp Dener isl[n_isl]=low+i; ++n_isl; 261c4b75bccSAlp Dener isa[n_isa]=low+i; ++n_isa; 262*76be6f4fSStefano Zampini } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - *bound_tol && g[i] < zero) { 263c4b75bccSAlp Dener /* Upper bounded variables */ 2642f75a4aaSAlp Dener isu[n_isu]=low+i; ++n_isu; 265c4b75bccSAlp Dener isa[n_isa]=low+i; ++n_isa; 266c4b75bccSAlp Dener } else { 267c4b75bccSAlp Dener /* Inactive variables */ 268c4b75bccSAlp Dener isi[n_isi]=low+i; ++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 } 3152f75a4aaSAlp Dener PetscFunctionReturn(0); 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 @*/ 336c4b75bccSAlp Dener PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, PetscReal scale, Vec S) 3372f75a4aaSAlp Dener { 3382f75a4aaSAlp Dener 3392f75a4aaSAlp Dener Vec step_lower, step_upper, step_fixed; 3402f75a4aaSAlp Dener Vec x_lower, x_upper; 3412f75a4aaSAlp Dener Vec bound_lower, bound_upper; 3422f75a4aaSAlp Dener 3432f75a4aaSAlp Dener PetscFunctionBegin; 3442f75a4aaSAlp Dener /* Adjust step for variables at the estimated lower bound */ 3452f75a4aaSAlp Dener if (active_lower) { 3469566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(S, active_lower, &step_lower)); 3479566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, active_lower, &x_lower)); 3489566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(XL, active_lower, &bound_lower)); 3499566063dSJacob Faibussowitsch PetscCall(VecCopy(bound_lower, step_lower)); 3509566063dSJacob Faibussowitsch PetscCall(VecAXPY(step_lower, -1.0, x_lower)); 3519566063dSJacob Faibussowitsch PetscCall(VecScale(step_lower, scale)); 3529566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(S, active_lower, &step_lower)); 3539566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, active_lower, &x_lower)); 3549566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(XL, active_lower, &bound_lower)); 3552f75a4aaSAlp Dener } 3562f75a4aaSAlp Dener 3572f75a4aaSAlp Dener /* Adjust step for the variables at the estimated upper bound */ 3582f75a4aaSAlp Dener if (active_upper) { 3599566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(S, active_upper, &step_upper)); 3609566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, active_upper, &x_upper)); 3619566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(XU, active_upper, &bound_upper)); 3629566063dSJacob Faibussowitsch PetscCall(VecCopy(bound_upper, step_upper)); 3639566063dSJacob Faibussowitsch PetscCall(VecAXPY(step_upper, -1.0, x_upper)); 3649566063dSJacob Faibussowitsch PetscCall(VecScale(step_upper, scale)); 3659566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(S, active_upper, &step_upper)); 3669566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, active_upper, &x_upper)); 3679566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(XU, active_upper, &bound_upper)); 3682f75a4aaSAlp Dener } 3692f75a4aaSAlp Dener 3702f75a4aaSAlp Dener /* Zero out step for fixed variables */ 3712f75a4aaSAlp Dener if (active_fixed) { 3729566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(S, active_fixed, &step_fixed)); 3739566063dSJacob Faibussowitsch PetscCall(VecSet(step_fixed, 0.0)); 3749566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(S, active_fixed, &step_fixed)); 3752f75a4aaSAlp Dener } 3762f75a4aaSAlp Dener PetscFunctionReturn(0); 3772f75a4aaSAlp Dener } 378c4b75bccSAlp Dener 379c4b75bccSAlp Dener /*@C 380b6a6cedcSAlp Dener TaoBoundSolution - Ensures that the solution vector is snapped into the bounds within a given tolerance. 381c4b75bccSAlp Dener 382f0fc11ceSJed Brown Collective on Vec 383f0fc11ceSJed Brown 384c4b75bccSAlp Dener Input Parameters: 385b6a6cedcSAlp Dener + X - solution vector 386b6a6cedcSAlp Dener . XL - lower bound vector 387c4b75bccSAlp Dener . XU - upper bound vector 388f0fc11ceSJed Brown - bound_tol - absolute tolerance in enforcing the bound 389c4b75bccSAlp Dener 390c4b75bccSAlp Dener Output Parameters: 391f0fc11ceSJed Brown + nDiff - total number of vector entries that have been bounded 392f0fc11ceSJed Brown - Xout - modified solution vector satisfying bounds to bound_tol 393b6a6cedcSAlp Dener 394b6a6cedcSAlp Dener Level: developer 395f0fc11ceSJed Brown 396db781477SPatrick Sanan .seealso: `TAOBNCG`, `TAOBNTL`, `TAOBNTR` 397c4b75bccSAlp Dener @*/ 3983b063059SAlp Dener PetscErrorCode TaoBoundSolution(Vec X, Vec XL, Vec XU, PetscReal bound_tol, PetscInt *nDiff, Vec Xout) 399c4b75bccSAlp Dener { 400c4b75bccSAlp Dener PetscInt i,n,low,high,nDiff_loc=0; 4013b063059SAlp Dener PetscScalar *xout; 4023b063059SAlp Dener const PetscScalar *x,*xl,*xu; 403c4b75bccSAlp Dener 404c4b75bccSAlp Dener PetscFunctionBegin; 405c4b75bccSAlp Dener PetscValidHeaderSpecific(X,VEC_CLASSID,1); 406*76be6f4fSStefano Zampini if (XL) PetscValidHeaderSpecific(XL,VEC_CLASSID,2); 407*76be6f4fSStefano Zampini if (XU) PetscValidHeaderSpecific(XU,VEC_CLASSID,3); 408064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(Xout,VEC_CLASSID,6); 409*76be6f4fSStefano Zampini if (!XL && !XU) { 410*76be6f4fSStefano Zampini PetscCall(VecCopy(X,Xout)); 411*76be6f4fSStefano Zampini *nDiff = 0.0; 412*76be6f4fSStefano Zampini PetscFunctionReturn(0); 413*76be6f4fSStefano Zampini } 414c4b75bccSAlp Dener PetscCheckSameType(X,1,XL,2); 415c4b75bccSAlp Dener PetscCheckSameType(X,1,XU,3); 416064a246eSJacob Faibussowitsch PetscCheckSameType(X,1,Xout,6); 417c4b75bccSAlp Dener PetscCheckSameComm(X,1,XL,2); 418c4b75bccSAlp Dener PetscCheckSameComm(X,1,XU,3); 419064a246eSJacob Faibussowitsch PetscCheckSameComm(X,1,Xout,6); 420c4b75bccSAlp Dener VecCheckSameSize(X,1,XL,2); 421c4b75bccSAlp Dener VecCheckSameSize(X,1,XU,3); 4223b063059SAlp Dener VecCheckSameSize(X,1,Xout,4); 423c4b75bccSAlp Dener 4249566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X,&low,&high)); 4259566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X,&n)); 426c4b75bccSAlp Dener if (n>0) { 4279566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 4289566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XL, &xl)); 4299566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XU, &xu)); 4309566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xout, &xout)); 431c4b75bccSAlp Dener 432c4b75bccSAlp Dener for (i=0;i<n;++i) { 433*76be6f4fSStefano Zampini if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + bound_tol) { 4343b063059SAlp Dener xout[i] = xl[i]; ++nDiff_loc; 435*76be6f4fSStefano Zampini } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - bound_tol) { 4363b063059SAlp Dener xout[i] = xu[i]; ++nDiff_loc; 437c4b75bccSAlp Dener } 438c4b75bccSAlp Dener } 439c4b75bccSAlp Dener 4409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 4419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XL, &xl)); 4429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XU, &xu)); 4439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xout, &xout)); 444c4b75bccSAlp Dener } 4451c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&nDiff_loc, nDiff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)X))); 446c4b75bccSAlp Dener PetscFunctionReturn(0); 447c4b75bccSAlp Dener } 448