xref: /petsc/src/tao/bound/utils/isutil.c (revision 1c2dc1cbabe5212f80cf309c4ca5a26f0cadc660)
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));
823c859ba3SBarry Smith       PetscCheck(nlocal <= (fhigh-flow),PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"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       }
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   PetscErrorCode ierr;
112a7e14dcfSSatish Balay   IS             iscomp;
11362675beeSAlp Dener   PetscBool      flg = PETSC_TRUE;
11453506e15SBarry Smith 
115a7e14dcfSSatish Balay   PetscFunctionBegin;
116a7e14dcfSSatish Balay   PetscValidHeaderSpecific(M,MAT_CLASSID,1);
117a7e14dcfSSatish Balay   PetscValidHeaderSpecific(is,IS_CLASSID,2);
1189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(Msub));
119a7e14dcfSSatish Balay   switch (subset_type) {
120a7e14dcfSSatish Balay   case TAO_SUBSET_SUBVEC:
1219566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub));
122a7e14dcfSSatish Balay     break;
123a7e14dcfSSatish Balay 
124a7e14dcfSSatish Balay   case TAO_SUBSET_MASK:
125a7e14dcfSSatish Balay     /* Get Reduced Hessian
126a7e14dcfSSatish Balay      Msub[i,j] = M[i,j] if i,j in Free_Local or i==j
127a7e14dcfSSatish Balay      Msub[i,j] = 0      if i!=j and i or j not in Free_Local
128a7e14dcfSSatish Balay      */
1299566063dSJacob Faibussowitsch     ierr = PetscObjectOptionsBegin((PetscObject)M);PetscCall(ierr);
1309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-overwrite_hessian","modify the existing hessian matrix when computing submatrices","TaoSubsetType",flg,&flg,NULL));
1319566063dSJacob Faibussowitsch     ierr = PetscOptionsEnd();PetscCall(ierr);
1328afaa268SBarry Smith     if (flg) {
1339566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(M, MAT_COPY_VALUES, Msub));
134a7e14dcfSSatish Balay     } else {
135a7e14dcfSSatish Balay       /* Act on hessian directly (default) */
1369566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)M));
137a7e14dcfSSatish Balay       *Msub = M;
138a7e14dcfSSatish Balay     }
139a7e14dcfSSatish Balay     /* Save the diagonal to temporary vector */
1409566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(*Msub,v1));
141a7e14dcfSSatish Balay 
142a7e14dcfSSatish Balay     /* Zero out rows and columns */
1439566063dSJacob Faibussowitsch     PetscCall(ISComplementVec(is,v1,&iscomp));
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay     /* Use v1 instead of 0 here because of PETSc bug */
1469566063dSJacob Faibussowitsch     PetscCall(MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1));
147a7e14dcfSSatish Balay 
1489566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iscomp));
149a7e14dcfSSatish Balay     break;
150a7e14dcfSSatish Balay   case TAO_SUBSET_MATRIXFREE:
1519566063dSJacob Faibussowitsch     PetscCall(ISComplementVec(is,v1,&iscomp));
1529566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrixFree(M,iscomp,iscomp,Msub));
1539566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iscomp));
154a7e14dcfSSatish Balay     break;
155a7e14dcfSSatish Balay   }
156a7e14dcfSSatish Balay   PetscFunctionReturn(0);
157a7e14dcfSSatish Balay }
1582f75a4aaSAlp Dener 
1592f75a4aaSAlp Dener /*@C
1602f75a4aaSAlp Dener   TaoEstimateActiveBounds - Generates index sets for variables at the lower and upper
1612f75a4aaSAlp Dener   bounds, as well as fixed variables where lower and upper bounds equal each other.
1622f75a4aaSAlp Dener 
1632f75a4aaSAlp Dener   Input Parameters:
1642f75a4aaSAlp Dener + X - solution vector
1652f75a4aaSAlp Dener . XL - lower bound vector
1662f75a4aaSAlp Dener . XU - upper bound vector
1672f75a4aaSAlp Dener . G - unprojected gradient
1680a4511e9SAlp Dener . S - step direction with which the active bounds will be estimated
169b6a6cedcSAlp Dener . W - work vector of type and size of X
1700a4511e9SAlp Dener - steplen - the step length at which the active bounds will be estimated (needs to be conservative)
1712f75a4aaSAlp Dener 
1722f75a4aaSAlp Dener   Output Parameters:
173b6a6cedcSAlp Dener + bound_tol - tolerance for for the bound estimation
1742f75a4aaSAlp Dener . active_lower - index set for active variables at the lower bound
1752f75a4aaSAlp Dener . active_upper - index set for active variables at the upper bound
1762f75a4aaSAlp Dener . active_fixed - index set for fixed variables
1772f75a4aaSAlp Dener . active - index set for all active variables
178b6a6cedcSAlp Dener - inactive - complementary index set for inactive variables
1790a4511e9SAlp Dener 
1800a4511e9SAlp Dener   Notes:
1810a4511e9SAlp Dener   This estimation is based on Bertsekas' method, with a built in diagonal scaling value of 1.0e-3.
1820a4511e9SAlp Dener 
183b6a6cedcSAlp Dener   Level: developer
1842f75a4aaSAlp Dener @*/
185c4b75bccSAlp Dener PetscErrorCode TaoEstimateActiveBounds(Vec X, Vec XL, Vec XU, Vec G, Vec S, Vec W, PetscReal steplen, PetscReal *bound_tol,
1862f75a4aaSAlp Dener                                        IS *active_lower, IS *active_upper, IS *active_fixed, IS *active, IS *inactive)
1872f75a4aaSAlp Dener {
1882f75a4aaSAlp Dener   PetscReal                    wnorm;
18989da521bSAlp Dener   PetscReal                    zero = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
190c4b75bccSAlp Dener   PetscInt                     i, n_isl=0, n_isu=0, n_isf=0, n_isa=0, n_isi=0;
191c4b75bccSAlp Dener   PetscInt                     N_isl, N_isu, N_isf, N_isa, N_isi;
192c4b75bccSAlp Dener   PetscInt                     n, low, high, nDiff;
193c4b75bccSAlp Dener   PetscInt                     *isl=NULL, *isu=NULL, *isf=NULL, *isa=NULL, *isi=NULL;
1942f75a4aaSAlp Dener   const PetscScalar            *xl, *xu, *x, *g;
1956519dc0cSAlp Dener   MPI_Comm                     comm = PetscObjectComm((PetscObject)X);
1962f75a4aaSAlp Dener 
1972f75a4aaSAlp Dener   PetscFunctionBegin;
198c4b75bccSAlp Dener   PetscValidHeaderSpecific(X,VEC_CLASSID,1);
199c4b75bccSAlp Dener   PetscValidHeaderSpecific(XL,VEC_CLASSID,2);
200c4b75bccSAlp Dener   PetscValidHeaderSpecific(XU,VEC_CLASSID,3);
201c4b75bccSAlp Dener   PetscValidHeaderSpecific(G,VEC_CLASSID,4);
202c4b75bccSAlp Dener   PetscValidHeaderSpecific(S,VEC_CLASSID,5);
203c4b75bccSAlp Dener   PetscValidHeaderSpecific(W,VEC_CLASSID,6);
204c4b75bccSAlp Dener 
205c4b75bccSAlp Dener   PetscValidType(X,1);
206c4b75bccSAlp Dener   PetscValidType(XL,2);
207c4b75bccSAlp Dener   PetscValidType(XU,3);
208c4b75bccSAlp Dener   PetscValidType(G,4);
209c4b75bccSAlp Dener   PetscValidType(S,5);
210c4b75bccSAlp Dener   PetscValidType(W,6);
211c4b75bccSAlp Dener   PetscCheckSameType(X,1,XL,2);
212c4b75bccSAlp Dener   PetscCheckSameType(X,1,XU,3);
213c4b75bccSAlp Dener   PetscCheckSameType(X,1,G,4);
214c4b75bccSAlp Dener   PetscCheckSameType(X,1,S,5);
215c4b75bccSAlp Dener   PetscCheckSameType(X,1,W,6);
216c4b75bccSAlp Dener   PetscCheckSameComm(X,1,XL,2);
217c4b75bccSAlp Dener   PetscCheckSameComm(X,1,XU,3);
218c4b75bccSAlp Dener   PetscCheckSameComm(X,1,G,4);
219c4b75bccSAlp Dener   PetscCheckSameComm(X,1,S,5);
220c4b75bccSAlp Dener   PetscCheckSameComm(X,1,W,6);
221c4b75bccSAlp Dener   VecCheckSameSize(X,1,XL,2);
222c4b75bccSAlp Dener   VecCheckSameSize(X,1,XU,3);
223c4b75bccSAlp Dener   VecCheckSameSize(X,1,G,4);
224c4b75bccSAlp Dener   VecCheckSameSize(X,1,S,5);
225c4b75bccSAlp Dener   VecCheckSameSize(X,1,W,6);
226c4b75bccSAlp Dener 
2272f75a4aaSAlp Dener   /* Update the tolerance for bound detection (this is based on Bertsekas' method) */
2289566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, W));
2299566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(W, steplen, 1.0, S));
2309566063dSJacob Faibussowitsch   PetscCall(TaoBoundSolution(W, XL, XU, 0.0, &nDiff, W));
2319566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(W, 1.0, -1.0, X));
2329566063dSJacob Faibussowitsch   PetscCall(VecNorm(W, NORM_2, &wnorm));
2332f75a4aaSAlp Dener   *bound_tol = PetscMin(*bound_tol, wnorm);
2342f75a4aaSAlp Dener 
2359566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &low, &high));
2369566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
2372f75a4aaSAlp Dener   if (n>0) {
2389566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, &x));
2399566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XL, &xl));
2409566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XU, &xu));
2419566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(G, &g));
2422f75a4aaSAlp Dener 
2432f75a4aaSAlp Dener     /* Loop over variables and categorize the indexes */
2449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isl));
2459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isu));
2469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isf));
2479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isa));
2489566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isi));
249c4b75bccSAlp Dener     for (i=0; i<n; ++i) {
25066b4d56fSTodd Munson       if (xl[i] == xu[i]) {
251c4b75bccSAlp Dener         /* Fixed variables */
2522f75a4aaSAlp Dener         isf[n_isf]=low+i; ++n_isf;
253c4b75bccSAlp Dener         isa[n_isa]=low+i; ++n_isa;
25489da521bSAlp Dener       } else if ((xl[i] > PETSC_NINFINITY) && (x[i] <= xl[i] + *bound_tol) && (g[i] > zero)) {
255c4b75bccSAlp Dener         /* Lower bounded variables */
2562f75a4aaSAlp Dener         isl[n_isl]=low+i; ++n_isl;
257c4b75bccSAlp Dener         isa[n_isa]=low+i; ++n_isa;
25889da521bSAlp Dener       } else if ((xu[i] < PETSC_INFINITY) && (x[i] >= xu[i] - *bound_tol) && (g[i] < zero)) {
259c4b75bccSAlp Dener         /* Upper bounded variables */
2602f75a4aaSAlp Dener         isu[n_isu]=low+i; ++n_isu;
261c4b75bccSAlp Dener         isa[n_isa]=low+i; ++n_isa;
262c4b75bccSAlp Dener       } else {
263c4b75bccSAlp Dener         /* Inactive variables */
264c4b75bccSAlp Dener         isi[n_isi]=low+i; ++n_isi;
2652f75a4aaSAlp Dener       }
2662f75a4aaSAlp Dener     }
2672f75a4aaSAlp Dener 
2689566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, &x));
2699566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XL, &xl));
2709566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XU, &xu));
2719566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(G, &g));
2722f75a4aaSAlp Dener   }
2732f75a4aaSAlp Dener 
274c4b75bccSAlp Dener   /* Clear all index sets */
2759566063dSJacob Faibussowitsch   PetscCall(ISDestroy(active_lower));
2769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(active_upper));
2779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(active_fixed));
2789566063dSJacob Faibussowitsch   PetscCall(ISDestroy(active));
2799566063dSJacob Faibussowitsch   PetscCall(ISDestroy(inactive));
280c4b75bccSAlp Dener 
281c4b75bccSAlp Dener   /* Collect global sizes */
282*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&n_isl, &N_isl, 1, MPIU_INT, MPI_SUM, comm));
283*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&n_isu, &N_isu, 1, MPIU_INT, MPI_SUM, comm));
284*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&n_isf, &N_isf, 1, MPIU_INT, MPI_SUM, comm));
285*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&n_isa, &N_isa, 1, MPIU_INT, MPI_SUM, comm));
286*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&n_isi, &N_isi, 1, MPIU_INT, MPI_SUM, comm));
287c4b75bccSAlp Dener 
288c4b75bccSAlp Dener   /* Create index set for lower bounded variables */
289c4b75bccSAlp Dener   if (N_isl > 0) {
2909566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isl, isl, PETSC_OWN_POINTER, active_lower));
2917e9273c4SAlp Dener   } else {
2929566063dSJacob Faibussowitsch     PetscCall(PetscFree(isl));
293c4b75bccSAlp Dener   }
294c4b75bccSAlp Dener   /* Create index set for upper bounded variables */
295c4b75bccSAlp Dener   if (N_isu > 0) {
2969566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isu, isu, PETSC_OWN_POINTER, active_upper));
2977e9273c4SAlp Dener   } else {
2989566063dSJacob Faibussowitsch     PetscCall(PetscFree(isu));
299c4b75bccSAlp Dener   }
300c4b75bccSAlp Dener   /* Create index set for fixed variables */
301c4b75bccSAlp Dener   if (N_isf > 0) {
3029566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isf, isf, PETSC_OWN_POINTER, active_fixed));
3037e9273c4SAlp Dener   } else {
3049566063dSJacob Faibussowitsch     PetscCall(PetscFree(isf));
305c4b75bccSAlp Dener   }
306c4b75bccSAlp Dener   /* Create index set for all actively bounded variables */
307c4b75bccSAlp Dener   if (N_isa > 0) {
3089566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isa, isa, PETSC_OWN_POINTER, active));
3097e9273c4SAlp Dener   } else {
3109566063dSJacob Faibussowitsch     PetscCall(PetscFree(isa));
311c4b75bccSAlp Dener   }
312c4b75bccSAlp Dener   /* Create index set for all inactive variables */
313c4b75bccSAlp Dener   if (N_isi > 0) {
3149566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isi, isi, PETSC_OWN_POINTER, inactive));
3157e9273c4SAlp Dener   } else {
3169566063dSJacob Faibussowitsch     PetscCall(PetscFree(isi));
317c4b75bccSAlp Dener   }
3182f75a4aaSAlp Dener 
3192f75a4aaSAlp Dener   /* Clean up and exit */
3202f75a4aaSAlp Dener   PetscFunctionReturn(0);
3212f75a4aaSAlp Dener }
3222f75a4aaSAlp Dener 
3232f75a4aaSAlp Dener /*@C
3242f75a4aaSAlp Dener   TaoBoundStep - Ensures the correct zero or adjusted step direction
3252f75a4aaSAlp Dener   values for active variables.
3262f75a4aaSAlp Dener 
3272f75a4aaSAlp Dener   Input Parameters:
3282f75a4aaSAlp Dener + X - solution vector
3292f75a4aaSAlp Dener . XL - lower bound vector
3302f75a4aaSAlp Dener . XU - upper bound vector
3312f75a4aaSAlp Dener . active_lower - index set for lower bounded active variables
3322f75a4aaSAlp Dener . active_upper - index set for lower bounded active variables
333b6a6cedcSAlp Dener . active_fixed - index set for fixed active variables
334b6a6cedcSAlp Dener - scale - amplification factor for the step that needs to be taken on actively bounded variables
3352f75a4aaSAlp Dener 
33697bb3fdcSJose E. Roman   Output Parameter:
3372f75a4aaSAlp Dener . S - step direction to be modified
338b6a6cedcSAlp Dener 
339b6a6cedcSAlp Dener   Level: developer
3402f75a4aaSAlp Dener @*/
341c4b75bccSAlp Dener PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, PetscReal scale, Vec S)
3422f75a4aaSAlp Dener {
3432f75a4aaSAlp Dener 
3442f75a4aaSAlp Dener   Vec                          step_lower, step_upper, step_fixed;
3452f75a4aaSAlp Dener   Vec                          x_lower, x_upper;
3462f75a4aaSAlp Dener   Vec                          bound_lower, bound_upper;
3472f75a4aaSAlp Dener 
3482f75a4aaSAlp Dener   PetscFunctionBegin;
3492f75a4aaSAlp Dener   /* Adjust step for variables at the estimated lower bound */
3502f75a4aaSAlp Dener   if (active_lower) {
3519566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(S, active_lower, &step_lower));
3529566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(X, active_lower, &x_lower));
3539566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(XL, active_lower, &bound_lower));
3549566063dSJacob Faibussowitsch     PetscCall(VecCopy(bound_lower, step_lower));
3559566063dSJacob Faibussowitsch     PetscCall(VecAXPY(step_lower, -1.0, x_lower));
3569566063dSJacob Faibussowitsch     PetscCall(VecScale(step_lower, scale));
3579566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(S, active_lower, &step_lower));
3589566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(X, active_lower, &x_lower));
3599566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(XL, active_lower, &bound_lower));
3602f75a4aaSAlp Dener   }
3612f75a4aaSAlp Dener 
3622f75a4aaSAlp Dener   /* Adjust step for the variables at the estimated upper bound */
3632f75a4aaSAlp Dener   if (active_upper) {
3649566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(S, active_upper, &step_upper));
3659566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(X, active_upper, &x_upper));
3669566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(XU, active_upper, &bound_upper));
3679566063dSJacob Faibussowitsch     PetscCall(VecCopy(bound_upper, step_upper));
3689566063dSJacob Faibussowitsch     PetscCall(VecAXPY(step_upper, -1.0, x_upper));
3699566063dSJacob Faibussowitsch     PetscCall(VecScale(step_upper, scale));
3709566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(S, active_upper, &step_upper));
3719566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(X, active_upper, &x_upper));
3729566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(XU, active_upper, &bound_upper));
3732f75a4aaSAlp Dener   }
3742f75a4aaSAlp Dener 
3752f75a4aaSAlp Dener   /* Zero out step for fixed variables */
3762f75a4aaSAlp Dener   if (active_fixed) {
3779566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(S, active_fixed, &step_fixed));
3789566063dSJacob Faibussowitsch     PetscCall(VecSet(step_fixed, 0.0));
3799566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(S, active_fixed, &step_fixed));
3802f75a4aaSAlp Dener   }
3812f75a4aaSAlp Dener   PetscFunctionReturn(0);
3822f75a4aaSAlp Dener }
383c4b75bccSAlp Dener 
384c4b75bccSAlp Dener /*@C
385b6a6cedcSAlp Dener   TaoBoundSolution - Ensures that the solution vector is snapped into the bounds within a given tolerance.
386c4b75bccSAlp Dener 
387f0fc11ceSJed Brown   Collective on Vec
388f0fc11ceSJed Brown 
389c4b75bccSAlp Dener   Input Parameters:
390b6a6cedcSAlp Dener + X - solution vector
391b6a6cedcSAlp Dener . XL - lower bound vector
392c4b75bccSAlp Dener . XU - upper bound vector
393f0fc11ceSJed Brown - bound_tol - absolute tolerance in enforcing the bound
394c4b75bccSAlp Dener 
395c4b75bccSAlp Dener   Output Parameters:
396f0fc11ceSJed Brown + nDiff - total number of vector entries that have been bounded
397f0fc11ceSJed Brown - Xout - modified solution vector satisfying bounds to bound_tol
398b6a6cedcSAlp Dener 
399b6a6cedcSAlp Dener   Level: developer
400f0fc11ceSJed Brown 
401f0fc11ceSJed Brown .seealso: TAOBNCG, TAOBNTL, TAOBNTR
402c4b75bccSAlp Dener @*/
4033b063059SAlp Dener PetscErrorCode TaoBoundSolution(Vec X, Vec XL, Vec XU, PetscReal bound_tol, PetscInt *nDiff, Vec Xout)
404c4b75bccSAlp Dener {
405c4b75bccSAlp Dener   PetscInt          i,n,low,high,nDiff_loc=0;
4063b063059SAlp Dener   PetscScalar       *xout;
4073b063059SAlp Dener   const PetscScalar *x,*xl,*xu;
408c4b75bccSAlp Dener 
409c4b75bccSAlp Dener   PetscFunctionBegin;
410c4b75bccSAlp Dener   PetscValidHeaderSpecific(X,VEC_CLASSID,1);
411c4b75bccSAlp Dener   PetscValidHeaderSpecific(XL,VEC_CLASSID,2);
412c4b75bccSAlp Dener   PetscValidHeaderSpecific(XU,VEC_CLASSID,3);
413064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Xout,VEC_CLASSID,6);
414c4b75bccSAlp Dener 
415c4b75bccSAlp Dener   PetscValidType(X,1);
416c4b75bccSAlp Dener   PetscValidType(XL,2);
417c4b75bccSAlp Dener   PetscValidType(XU,3);
418064a246eSJacob Faibussowitsch   PetscValidType(Xout,6);
419c4b75bccSAlp Dener   PetscCheckSameType(X,1,XL,2);
420c4b75bccSAlp Dener   PetscCheckSameType(X,1,XU,3);
421064a246eSJacob Faibussowitsch   PetscCheckSameType(X,1,Xout,6);
422c4b75bccSAlp Dener   PetscCheckSameComm(X,1,XL,2);
423c4b75bccSAlp Dener   PetscCheckSameComm(X,1,XU,3);
424064a246eSJacob Faibussowitsch   PetscCheckSameComm(X,1,Xout,6);
425c4b75bccSAlp Dener   VecCheckSameSize(X,1,XL,2);
426c4b75bccSAlp Dener   VecCheckSameSize(X,1,XU,3);
4273b063059SAlp Dener   VecCheckSameSize(X,1,Xout,4);
428c4b75bccSAlp Dener 
4299566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X,&low,&high));
4309566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&n));
431c4b75bccSAlp Dener   if (n>0) {
4329566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, &x));
4339566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XL, &xl));
4349566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XU, &xu));
4359566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Xout, &xout));
436c4b75bccSAlp Dener 
437c4b75bccSAlp Dener     for (i=0;i<n;++i) {
43889da521bSAlp Dener       if ((xl[i] > PETSC_NINFINITY) && (x[i] <= xl[i] + bound_tol)) {
4393b063059SAlp Dener         xout[i] = xl[i]; ++nDiff_loc;
44089da521bSAlp Dener       } else if ((xu[i] < PETSC_INFINITY) && (x[i] >= xu[i] - bound_tol)) {
4413b063059SAlp Dener         xout[i] = xu[i]; ++nDiff_loc;
442c4b75bccSAlp Dener       }
443c4b75bccSAlp Dener     }
444c4b75bccSAlp Dener 
4459566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, &x));
4469566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XL, &xl));
4479566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XU, &xu));
4489566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Xout, &xout));
449c4b75bccSAlp Dener   }
450*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&nDiff_loc, nDiff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)X)));
451c4b75bccSAlp Dener   PetscFunctionReturn(0);
452c4b75bccSAlp Dener }
453