xref: /petsc/src/tao/bound/utils/isutil.c (revision f0fc11cebb1bb284829732915f9e84cabc170c2f)
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 
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   PetscErrorCode ierr;
26a7e14dcfSSatish Balay   PetscInt       nfull,nreduced,nreduced_local,rlow,rhigh,flow,fhigh;
27a7e14dcfSSatish Balay   PetscInt       i,nlocal;
28a7e14dcfSSatish Balay   PetscReal      *fv,*rv;
29a7e14dcfSSatish Balay   const PetscInt *s;
30a7e14dcfSSatish Balay   IS             ident;
31a7e14dcfSSatish Balay   VecType        vtype;
32a7e14dcfSSatish Balay   VecScatter     scatter;
33a7e14dcfSSatish Balay   MPI_Comm       comm;
34a7e14dcfSSatish Balay 
35a7e14dcfSSatish Balay   PetscFunctionBegin;
36a7e14dcfSSatish Balay   PetscValidHeaderSpecific(vfull,VEC_CLASSID,1);
37a7e14dcfSSatish Balay   PetscValidHeaderSpecific(is,IS_CLASSID,2);
38a7e14dcfSSatish Balay 
39a7e14dcfSSatish Balay   ierr = VecGetSize(vfull, &nfull);CHKERRQ(ierr);
40a7e14dcfSSatish Balay   ierr = ISGetSize(is, &nreduced);CHKERRQ(ierr);
41a7e14dcfSSatish Balay 
42a7e14dcfSSatish Balay   if (nreduced == nfull) {
43a7e14dcfSSatish Balay     ierr = VecDestroy(vreduced);CHKERRQ(ierr);
44a7e14dcfSSatish Balay     ierr = VecDuplicate(vfull,vreduced);CHKERRQ(ierr);
45a7e14dcfSSatish Balay     ierr = VecCopy(vfull,*vreduced);CHKERRQ(ierr);
46a7e14dcfSSatish Balay   } else {
47a7e14dcfSSatish Balay     switch (reduced_type) {
48a7e14dcfSSatish Balay     case TAO_SUBSET_SUBVEC:
49a7e14dcfSSatish Balay       ierr = VecGetType(vfull,&vtype);CHKERRQ(ierr);
50a7e14dcfSSatish Balay       ierr = VecGetOwnershipRange(vfull,&flow,&fhigh);CHKERRQ(ierr);
51a7e14dcfSSatish Balay       ierr = ISGetLocalSize(is,&nreduced_local);CHKERRQ(ierr);
52a7e14dcfSSatish Balay       ierr = PetscObjectGetComm((PetscObject)vfull,&comm);CHKERRQ(ierr);
53a7e14dcfSSatish Balay       if (*vreduced) {
54a7e14dcfSSatish Balay         ierr = VecDestroy(vreduced);CHKERRQ(ierr);
55a7e14dcfSSatish Balay       }
56a7e14dcfSSatish Balay       ierr = VecCreate(comm,vreduced);CHKERRQ(ierr);
57a7e14dcfSSatish Balay       ierr = VecSetType(*vreduced,vtype);CHKERRQ(ierr);
58a7e14dcfSSatish Balay 
59a7e14dcfSSatish Balay       ierr = VecSetSizes(*vreduced,nreduced_local,nreduced);CHKERRQ(ierr);
60a7e14dcfSSatish Balay       ierr = VecGetOwnershipRange(*vreduced,&rlow,&rhigh);CHKERRQ(ierr);
61a7e14dcfSSatish Balay       ierr = ISCreateStride(comm,nreduced_local,rlow,1,&ident);CHKERRQ(ierr);
629448b7f1SJunchao Zhang       ierr = VecScatterCreate(vfull,is,*vreduced,ident,&scatter);CHKERRQ(ierr);
63a7e14dcfSSatish Balay       ierr = VecScatterBegin(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
64a7e14dcfSSatish Balay       ierr = VecScatterEnd(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
65a7e14dcfSSatish Balay       ierr = VecScatterDestroy(&scatter);CHKERRQ(ierr);
66a7e14dcfSSatish Balay       ierr = ISDestroy(&ident);CHKERRQ(ierr);
67a7e14dcfSSatish Balay       break;
68a7e14dcfSSatish Balay 
69a7e14dcfSSatish Balay     case TAO_SUBSET_MASK:
70a7e14dcfSSatish Balay     case TAO_SUBSET_MATRIXFREE:
71a7e14dcfSSatish Balay       /* vr[i] = vf[i]   if i in is
72a7e14dcfSSatish Balay        vr[i] = 0       otherwise */
736c4ed002SBarry Smith       if (!*vreduced) {
74a7e14dcfSSatish Balay         ierr = VecDuplicate(vfull,vreduced);CHKERRQ(ierr);
75a7e14dcfSSatish Balay       }
76a7e14dcfSSatish Balay 
77a7e14dcfSSatish Balay       ierr = VecSet(*vreduced,maskvalue);CHKERRQ(ierr);
78a7e14dcfSSatish Balay       ierr = ISGetLocalSize(is,&nlocal);CHKERRQ(ierr);
79a7e14dcfSSatish Balay       ierr = VecGetOwnershipRange(vfull,&flow,&fhigh);CHKERRQ(ierr);
80a7e14dcfSSatish Balay       ierr = VecGetArray(vfull,&fv);CHKERRQ(ierr);
81a7e14dcfSSatish Balay       ierr = VecGetArray(*vreduced,&rv);CHKERRQ(ierr);
82a7e14dcfSSatish Balay       ierr = ISGetIndices(is,&s);CHKERRQ(ierr);
8353506e15SBarry Smith       if (nlocal > (fhigh-flow)) SETERRQ2(PETSC_COMM_WORLD,1,"IS local size %d > Vec local size %d",nlocal,fhigh-flow);
84c4b75bccSAlp Dener       for (i=0;i<nlocal;++i) {
85a7e14dcfSSatish Balay         rv[s[i]-flow] = fv[s[i]-flow];
86a7e14dcfSSatish Balay       }
87a7e14dcfSSatish Balay       ierr = ISRestoreIndices(is,&s);CHKERRQ(ierr);
88a7e14dcfSSatish Balay       ierr = VecRestoreArray(vfull,&fv);CHKERRQ(ierr);
89a7e14dcfSSatish Balay       ierr = VecRestoreArray(*vreduced,&rv);CHKERRQ(ierr);
90a7e14dcfSSatish Balay       break;
91a7e14dcfSSatish Balay     }
92a7e14dcfSSatish Balay   }
93a7e14dcfSSatish Balay   PetscFunctionReturn(0);
94a7e14dcfSSatish Balay }
95a7e14dcfSSatish Balay 
96b98f30f2SJason Sarich /*@C
97b98f30f2SJason Sarich   TaoMatGetSubMat - Gets a submatrix using the IS
98a7e14dcfSSatish Balay 
99a7e14dcfSSatish Balay   Input Parameters:
100a7e14dcfSSatish Balay + M - the full matrix (n x n)
101a7e14dcfSSatish Balay . is - the index set for the submatrix (both row and column index sets need to be the same)
102a7e14dcfSSatish Balay . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option
103b6a6cedcSAlp Dener - subset_type <TAO_SUBSET_SUBVEC,TAO_SUBSET_MASK,TAO_SUBSET_MATRIXFREE> - the method TAO is using for subsetting
104a7e14dcfSSatish Balay 
105a7e14dcfSSatish Balay   Output Parameters:
106a7e14dcfSSatish Balay . Msub - the submatrix
107b6a6cedcSAlp Dener 
108b6a6cedcSAlp Dener   Level: developer
109a7e14dcfSSatish Balay @*/
110b98f30f2SJason Sarich PetscErrorCode TaoMatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub)
111a7e14dcfSSatish Balay {
112a7e14dcfSSatish Balay   PetscErrorCode ierr;
113a7e14dcfSSatish Balay   IS             iscomp;
11462675beeSAlp Dener   PetscBool      flg = PETSC_TRUE;
11553506e15SBarry Smith 
116a7e14dcfSSatish Balay   PetscFunctionBegin;
117a7e14dcfSSatish Balay   PetscValidHeaderSpecific(M,MAT_CLASSID,1);
118a7e14dcfSSatish Balay   PetscValidHeaderSpecific(is,IS_CLASSID,2);
119a7e14dcfSSatish Balay   ierr = MatDestroy(Msub);CHKERRQ(ierr);
120a7e14dcfSSatish Balay   switch (subset_type) {
121a7e14dcfSSatish Balay   case TAO_SUBSET_SUBVEC:
1227dae84e0SHong Zhang     ierr = MatCreateSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub);CHKERRQ(ierr);
123a7e14dcfSSatish Balay     break;
124a7e14dcfSSatish Balay 
125a7e14dcfSSatish Balay   case TAO_SUBSET_MASK:
126a7e14dcfSSatish Balay     /* Get Reduced Hessian
127a7e14dcfSSatish Balay      Msub[i,j] = M[i,j] if i,j in Free_Local or i==j
128a7e14dcfSSatish Balay      Msub[i,j] = 0      if i!=j and i or j not in Free_Local
129a7e14dcfSSatish Balay      */
1301a1499c8SBarry Smith     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)M),NULL,NULL,NULL);CHKERRQ(ierr);
13162675beeSAlp Dener     ierr = PetscOptionsBool("-overwrite_hessian","modify the existing hessian matrix when computing submatrices","TaoSubsetType",flg,&flg,NULL);CHKERRQ(ierr);
132302440fdSBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1338afaa268SBarry Smith     if (flg) {
134a7e14dcfSSatish Balay       ierr = MatDuplicate(M, MAT_COPY_VALUES, Msub);CHKERRQ(ierr);
135a7e14dcfSSatish Balay     } else {
136a7e14dcfSSatish Balay       /* Act on hessian directly (default) */
137a7e14dcfSSatish Balay       ierr = PetscObjectReference((PetscObject)M);CHKERRQ(ierr);
138a7e14dcfSSatish Balay       *Msub = M;
139a7e14dcfSSatish Balay     }
140a7e14dcfSSatish Balay     /* Save the diagonal to temporary vector */
141a7e14dcfSSatish Balay     ierr = MatGetDiagonal(*Msub,v1);CHKERRQ(ierr);
142a7e14dcfSSatish Balay 
143a7e14dcfSSatish Balay     /* Zero out rows and columns */
1444473680cSBarry Smith     ierr = ISComplementVec(is,v1,&iscomp);CHKERRQ(ierr);
145a7e14dcfSSatish Balay 
146a7e14dcfSSatish Balay     /* Use v1 instead of 0 here because of PETSc bug */
147a7e14dcfSSatish Balay     ierr = MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1);CHKERRQ(ierr);
148a7e14dcfSSatish Balay 
149a7e14dcfSSatish Balay     ierr = ISDestroy(&iscomp);CHKERRQ(ierr);
150a7e14dcfSSatish Balay     break;
151a7e14dcfSSatish Balay   case TAO_SUBSET_MATRIXFREE:
1524473680cSBarry Smith     ierr = ISComplementVec(is,v1,&iscomp);CHKERRQ(ierr);
153a7e14dcfSSatish Balay     ierr = MatCreateSubMatrixFree(M,iscomp,iscomp,Msub);CHKERRQ(ierr);
154a7e14dcfSSatish Balay     ierr = ISDestroy(&iscomp);CHKERRQ(ierr);
155a7e14dcfSSatish Balay     break;
156a7e14dcfSSatish Balay   }
157a7e14dcfSSatish Balay   PetscFunctionReturn(0);
158a7e14dcfSSatish Balay }
1592f75a4aaSAlp Dener 
1602f75a4aaSAlp Dener /*@C
1612f75a4aaSAlp Dener   TaoEstimateActiveBounds - Generates index sets for variables at the lower and upper
1622f75a4aaSAlp Dener   bounds, as well as fixed variables where lower and upper bounds equal each other.
1632f75a4aaSAlp Dener 
1642f75a4aaSAlp Dener   Input Parameters:
1652f75a4aaSAlp Dener + X - solution vector
1662f75a4aaSAlp Dener . XL - lower bound vector
1672f75a4aaSAlp Dener . XU - upper bound vector
1682f75a4aaSAlp Dener . G - unprojected gradient
1690a4511e9SAlp Dener . S - step direction with which the active bounds will be estimated
170b6a6cedcSAlp Dener . W - work vector of type and size of X
1710a4511e9SAlp Dener - steplen - the step length at which the active bounds will be estimated (needs to be conservative)
1722f75a4aaSAlp Dener 
1732f75a4aaSAlp Dener   Output Parameters:
174b6a6cedcSAlp Dener + bound_tol - tolerance for for the bound estimation
1752f75a4aaSAlp Dener . active_lower - index set for active variables at the lower bound
1762f75a4aaSAlp Dener . active_upper - index set for active variables at the upper bound
1772f75a4aaSAlp Dener . active_fixed - index set for fixed variables
1782f75a4aaSAlp Dener . active - index set for all active variables
179b6a6cedcSAlp Dener - inactive - complementary index set for inactive variables
1800a4511e9SAlp Dener 
1810a4511e9SAlp Dener   Notes:
1820a4511e9SAlp Dener   This estimation is based on Bertsekas' method, with a built in diagonal scaling value of 1.0e-3.
1830a4511e9SAlp Dener 
184b6a6cedcSAlp Dener   Level: developer
1852f75a4aaSAlp Dener @*/
186c4b75bccSAlp Dener PetscErrorCode TaoEstimateActiveBounds(Vec X, Vec XL, Vec XU, Vec G, Vec S, Vec W, PetscReal steplen, PetscReal *bound_tol,
1872f75a4aaSAlp Dener                                        IS *active_lower, IS *active_upper, IS *active_fixed, IS *active, IS *inactive)
1882f75a4aaSAlp Dener {
1892f75a4aaSAlp Dener   PetscErrorCode               ierr;
1902f75a4aaSAlp Dener   PetscReal                    wnorm;
19189da521bSAlp Dener   PetscReal                    zero = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
192c4b75bccSAlp Dener   PetscInt                     i, n_isl=0, n_isu=0, n_isf=0, n_isa=0, n_isi=0;
193c4b75bccSAlp Dener   PetscInt                     N_isl, N_isu, N_isf, N_isa, N_isi;
194c4b75bccSAlp Dener   PetscInt                     n, low, high, nDiff;
195c4b75bccSAlp Dener   PetscInt                     *isl=NULL, *isu=NULL, *isf=NULL, *isa=NULL, *isi=NULL;
1962f75a4aaSAlp Dener   const PetscScalar            *xl, *xu, *x, *g;
1976519dc0cSAlp Dener   MPI_Comm                     comm = PetscObjectComm((PetscObject)X);
1982f75a4aaSAlp Dener 
1992f75a4aaSAlp Dener   PetscFunctionBegin;
200c4b75bccSAlp Dener   PetscValidHeaderSpecific(X,VEC_CLASSID,1);
201c4b75bccSAlp Dener   PetscValidHeaderSpecific(XL,VEC_CLASSID,2);
202c4b75bccSAlp Dener   PetscValidHeaderSpecific(XU,VEC_CLASSID,3);
203c4b75bccSAlp Dener   PetscValidHeaderSpecific(G,VEC_CLASSID,4);
204c4b75bccSAlp Dener   PetscValidHeaderSpecific(S,VEC_CLASSID,5);
205c4b75bccSAlp Dener   PetscValidHeaderSpecific(W,VEC_CLASSID,6);
206c4b75bccSAlp Dener 
207c4b75bccSAlp Dener   PetscValidType(X,1);
208c4b75bccSAlp Dener   PetscValidType(XL,2);
209c4b75bccSAlp Dener   PetscValidType(XU,3);
210c4b75bccSAlp Dener   PetscValidType(G,4);
211c4b75bccSAlp Dener   PetscValidType(S,5);
212c4b75bccSAlp Dener   PetscValidType(W,6);
213c4b75bccSAlp Dener   PetscCheckSameType(X,1,XL,2);
214c4b75bccSAlp Dener   PetscCheckSameType(X,1,XU,3);
215c4b75bccSAlp Dener   PetscCheckSameType(X,1,G,4);
216c4b75bccSAlp Dener   PetscCheckSameType(X,1,S,5);
217c4b75bccSAlp Dener   PetscCheckSameType(X,1,W,6);
218c4b75bccSAlp Dener   PetscCheckSameComm(X,1,XL,2);
219c4b75bccSAlp Dener   PetscCheckSameComm(X,1,XU,3);
220c4b75bccSAlp Dener   PetscCheckSameComm(X,1,G,4);
221c4b75bccSAlp Dener   PetscCheckSameComm(X,1,S,5);
222c4b75bccSAlp Dener   PetscCheckSameComm(X,1,W,6);
223c4b75bccSAlp Dener   VecCheckSameSize(X,1,XL,2);
224c4b75bccSAlp Dener   VecCheckSameSize(X,1,XU,3);
225c4b75bccSAlp Dener   VecCheckSameSize(X,1,G,4);
226c4b75bccSAlp Dener   VecCheckSameSize(X,1,S,5);
227c4b75bccSAlp Dener   VecCheckSameSize(X,1,W,6);
228c4b75bccSAlp Dener 
2292f75a4aaSAlp Dener   /* Update the tolerance for bound detection (this is based on Bertsekas' method) */
230c4b75bccSAlp Dener   ierr = VecCopy(X, W);CHKERRQ(ierr);
231c4b75bccSAlp Dener   ierr = VecAXPBY(W, steplen, 1.0, S);CHKERRQ(ierr);
2323b063059SAlp Dener   ierr = TaoBoundSolution(W, XL, XU, 0.0, &nDiff, W);CHKERRQ(ierr);
2332f75a4aaSAlp Dener   ierr = VecAXPBY(W, 1.0, -1.0, X);CHKERRQ(ierr);
2342f75a4aaSAlp Dener   ierr = VecNorm(W, NORM_2, &wnorm);CHKERRQ(ierr);
2352f75a4aaSAlp Dener   *bound_tol = PetscMin(*bound_tol, wnorm);
2362f75a4aaSAlp Dener 
2372f75a4aaSAlp Dener   ierr = VecGetOwnershipRange(X, &low, &high);CHKERRQ(ierr);
2382f75a4aaSAlp Dener   ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr);
2392f75a4aaSAlp Dener   if (n>0){
2402f75a4aaSAlp Dener     ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
2412f75a4aaSAlp Dener     ierr = VecGetArrayRead(XL, &xl);CHKERRQ(ierr);
2422f75a4aaSAlp Dener     ierr = VecGetArrayRead(XU, &xu);CHKERRQ(ierr);
2432f75a4aaSAlp Dener     ierr = VecGetArrayRead(G, &g);CHKERRQ(ierr);
2442f75a4aaSAlp Dener 
2452f75a4aaSAlp Dener     /* Loop over variables and categorize the indexes */
2462f75a4aaSAlp Dener     ierr = PetscMalloc1(n, &isl);CHKERRQ(ierr);
2472f75a4aaSAlp Dener     ierr = PetscMalloc1(n, &isu);CHKERRQ(ierr);
2482f75a4aaSAlp Dener     ierr = PetscMalloc1(n, &isf);CHKERRQ(ierr);
249c4b75bccSAlp Dener     ierr = PetscMalloc1(n, &isa);CHKERRQ(ierr);
250c4b75bccSAlp Dener     ierr = PetscMalloc1(n, &isi);CHKERRQ(ierr);
251c4b75bccSAlp Dener     for (i=0; i<n; ++i) {
25266b4d56fSTodd Munson       if (xl[i] == xu[i]) {
253c4b75bccSAlp Dener         /* Fixed variables */
2542f75a4aaSAlp Dener         isf[n_isf]=low+i; ++n_isf;
255c4b75bccSAlp Dener         isa[n_isa]=low+i; ++n_isa;
25689da521bSAlp Dener       } else if ((xl[i] > PETSC_NINFINITY) && (x[i] <= xl[i] + *bound_tol) && (g[i] > zero)) {
257c4b75bccSAlp Dener         /* Lower bounded variables */
2582f75a4aaSAlp Dener         isl[n_isl]=low+i; ++n_isl;
259c4b75bccSAlp Dener         isa[n_isa]=low+i; ++n_isa;
26089da521bSAlp Dener       } else if ((xu[i] < PETSC_INFINITY) && (x[i] >= xu[i] - *bound_tol) && (g[i] < zero)) {
261c4b75bccSAlp Dener         /* Upper bounded variables */
2622f75a4aaSAlp Dener         isu[n_isu]=low+i; ++n_isu;
263c4b75bccSAlp Dener         isa[n_isa]=low+i; ++n_isa;
264c4b75bccSAlp Dener       } else {
265c4b75bccSAlp Dener         /* Inactive variables */
266c4b75bccSAlp Dener         isi[n_isi]=low+i; ++n_isi;
2672f75a4aaSAlp Dener       }
2682f75a4aaSAlp Dener     }
2692f75a4aaSAlp Dener 
2702f75a4aaSAlp Dener     ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);
2712f75a4aaSAlp Dener     ierr = VecRestoreArrayRead(XL, &xl);CHKERRQ(ierr);
2722f75a4aaSAlp Dener     ierr = VecRestoreArrayRead(XU, &xu);CHKERRQ(ierr);
2732f75a4aaSAlp Dener     ierr = VecRestoreArrayRead(G, &g);CHKERRQ(ierr);
2742f75a4aaSAlp Dener   }
2752f75a4aaSAlp Dener 
276c4b75bccSAlp Dener   /* Clear all index sets */
2772f75a4aaSAlp Dener   ierr = ISDestroy(active_lower);CHKERRQ(ierr);
2782f75a4aaSAlp Dener   ierr = ISDestroy(active_upper);CHKERRQ(ierr);
2792f75a4aaSAlp Dener   ierr = ISDestroy(active_fixed);CHKERRQ(ierr);
2802f75a4aaSAlp Dener   ierr = ISDestroy(active);CHKERRQ(ierr);
2812f75a4aaSAlp Dener   ierr = ISDestroy(inactive);CHKERRQ(ierr);
282c4b75bccSAlp Dener 
283c4b75bccSAlp Dener   /* Collect global sizes */
2846519dc0cSAlp Dener   ierr = MPIU_Allreduce(&n_isl, &N_isl, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
2856519dc0cSAlp Dener   ierr = MPIU_Allreduce(&n_isu, &N_isu, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
2866519dc0cSAlp Dener   ierr = MPIU_Allreduce(&n_isf, &N_isf, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
2876519dc0cSAlp Dener   ierr = MPIU_Allreduce(&n_isa, &N_isa, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
2886519dc0cSAlp Dener   ierr = MPIU_Allreduce(&n_isi, &N_isi, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
289c4b75bccSAlp Dener 
290c4b75bccSAlp Dener   /* Create index set for lower bounded variables */
291c4b75bccSAlp Dener   if (N_isl > 0) {
2926519dc0cSAlp Dener     ierr = ISCreateGeneral(comm, n_isl, isl, PETSC_OWN_POINTER, active_lower);CHKERRQ(ierr);
2937e9273c4SAlp Dener   } else {
2947e9273c4SAlp Dener     ierr = PetscFree(isl);CHKERRQ(ierr);
295c4b75bccSAlp Dener   }
296c4b75bccSAlp Dener   /* Create index set for upper bounded variables */
297c4b75bccSAlp Dener   if (N_isu > 0) {
2986519dc0cSAlp Dener     ierr = ISCreateGeneral(comm, n_isu, isu, PETSC_OWN_POINTER, active_upper);CHKERRQ(ierr);
2997e9273c4SAlp Dener   } else {
3007e9273c4SAlp Dener     ierr = PetscFree(isu);CHKERRQ(ierr);
301c4b75bccSAlp Dener   }
302c4b75bccSAlp Dener   /* Create index set for fixed variables */
303c4b75bccSAlp Dener   if (N_isf > 0) {
3046519dc0cSAlp Dener     ierr = ISCreateGeneral(comm, n_isf, isf, PETSC_OWN_POINTER, active_fixed);CHKERRQ(ierr);
3057e9273c4SAlp Dener   } else {
3067e9273c4SAlp Dener     ierr = PetscFree(isf);CHKERRQ(ierr);
307c4b75bccSAlp Dener   }
308c4b75bccSAlp Dener   /* Create index set for all actively bounded variables */
309c4b75bccSAlp Dener   if (N_isa > 0) {
3106519dc0cSAlp Dener     ierr = ISCreateGeneral(comm, n_isa, isa, PETSC_OWN_POINTER, active);CHKERRQ(ierr);
3117e9273c4SAlp Dener   } else {
3127e9273c4SAlp Dener     ierr = PetscFree(isa);CHKERRQ(ierr);
313c4b75bccSAlp Dener   }
314c4b75bccSAlp Dener   /* Create index set for all inactive variables */
315c4b75bccSAlp Dener   if (N_isi > 0) {
3166519dc0cSAlp Dener     ierr = ISCreateGeneral(comm, n_isi, isi, PETSC_OWN_POINTER, inactive);CHKERRQ(ierr);
3177e9273c4SAlp Dener   } else {
3187e9273c4SAlp Dener     ierr = PetscFree(isi);CHKERRQ(ierr);
319c4b75bccSAlp Dener   }
3202f75a4aaSAlp Dener 
3212f75a4aaSAlp Dener   /* Clean up and exit */
3222f75a4aaSAlp Dener   PetscFunctionReturn(0);
3232f75a4aaSAlp Dener }
3242f75a4aaSAlp Dener 
3252f75a4aaSAlp Dener /*@C
3262f75a4aaSAlp Dener   TaoBoundStep - Ensures the correct zero or adjusted step direction
3272f75a4aaSAlp Dener   values for active variables.
3282f75a4aaSAlp Dener 
3292f75a4aaSAlp Dener   Input Parameters:
3302f75a4aaSAlp Dener + X - solution vector
3312f75a4aaSAlp Dener . XL - lower bound vector
3322f75a4aaSAlp Dener . XU - upper bound vector
3332f75a4aaSAlp Dener . active_lower - index set for lower bounded active variables
3342f75a4aaSAlp Dener . active_upper - index set for lower bounded active variables
335b6a6cedcSAlp Dener . active_fixed - index set for fixed active variables
336b6a6cedcSAlp Dener - scale - amplification factor for the step that needs to be taken on actively bounded variables
3372f75a4aaSAlp Dener 
3382f75a4aaSAlp Dener   Output Parameters:
3392f75a4aaSAlp Dener . S - step direction to be modified
340b6a6cedcSAlp Dener 
341b6a6cedcSAlp Dener   Level: developer
3422f75a4aaSAlp Dener @*/
343c4b75bccSAlp Dener PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, PetscReal scale, Vec S)
3442f75a4aaSAlp Dener {
3452f75a4aaSAlp Dener   PetscErrorCode               ierr;
3462f75a4aaSAlp Dener 
3472f75a4aaSAlp Dener   Vec                          step_lower, step_upper, step_fixed;
3482f75a4aaSAlp Dener   Vec                          x_lower, x_upper;
3492f75a4aaSAlp Dener   Vec                          bound_lower, bound_upper;
3502f75a4aaSAlp Dener 
3512f75a4aaSAlp Dener   PetscFunctionBegin;
3522f75a4aaSAlp Dener   /* Adjust step for variables at the estimated lower bound */
3532f75a4aaSAlp Dener   if (active_lower) {
3542f75a4aaSAlp Dener     ierr = VecGetSubVector(S, active_lower, &step_lower);CHKERRQ(ierr);
3552f75a4aaSAlp Dener     ierr = VecGetSubVector(X, active_lower, &x_lower);CHKERRQ(ierr);
3562f75a4aaSAlp Dener     ierr = VecGetSubVector(XL, active_lower, &bound_lower);CHKERRQ(ierr);
3572f75a4aaSAlp Dener     ierr = VecCopy(bound_lower, step_lower);CHKERRQ(ierr);
3582f75a4aaSAlp Dener     ierr = VecAXPY(step_lower, -1.0, x_lower);CHKERRQ(ierr);
359c4b75bccSAlp Dener     ierr = VecScale(step_lower, scale);CHKERRQ(ierr);
3602f75a4aaSAlp Dener     ierr = VecRestoreSubVector(S, active_lower, &step_lower);CHKERRQ(ierr);
3612f75a4aaSAlp Dener     ierr = VecRestoreSubVector(X, active_lower, &x_lower);CHKERRQ(ierr);
3622f75a4aaSAlp Dener     ierr = VecRestoreSubVector(XL, active_lower, &bound_lower);CHKERRQ(ierr);
3632f75a4aaSAlp Dener   }
3642f75a4aaSAlp Dener 
3652f75a4aaSAlp Dener   /* Adjust step for the variables at the estimated upper bound */
3662f75a4aaSAlp Dener   if (active_upper) {
3672f75a4aaSAlp Dener     ierr = VecGetSubVector(S, active_upper, &step_upper);CHKERRQ(ierr);
3682f75a4aaSAlp Dener     ierr = VecGetSubVector(X, active_upper, &x_upper);CHKERRQ(ierr);
3692f75a4aaSAlp Dener     ierr = VecGetSubVector(XU, active_upper, &bound_upper);CHKERRQ(ierr);
3702f75a4aaSAlp Dener     ierr = VecCopy(bound_upper, step_upper);CHKERRQ(ierr);
3712f75a4aaSAlp Dener     ierr = VecAXPY(step_upper, -1.0, x_upper);CHKERRQ(ierr);
372c4b75bccSAlp Dener     ierr = VecScale(step_upper, scale);CHKERRQ(ierr);
3732f75a4aaSAlp Dener     ierr = VecRestoreSubVector(S, active_upper, &step_upper);CHKERRQ(ierr);
3742f75a4aaSAlp Dener     ierr = VecRestoreSubVector(X, active_upper, &x_upper);CHKERRQ(ierr);
3752f75a4aaSAlp Dener     ierr = VecRestoreSubVector(XU, active_upper, &bound_upper);CHKERRQ(ierr);
3762f75a4aaSAlp Dener   }
3772f75a4aaSAlp Dener 
3782f75a4aaSAlp Dener   /* Zero out step for fixed variables */
3792f75a4aaSAlp Dener   if (active_fixed) {
3802f75a4aaSAlp Dener     ierr = VecGetSubVector(S, active_fixed, &step_fixed);CHKERRQ(ierr);
3812f75a4aaSAlp Dener     ierr = VecSet(step_fixed, 0.0);CHKERRQ(ierr);
3822f75a4aaSAlp Dener     ierr = VecRestoreSubVector(S, active_fixed, &step_fixed);CHKERRQ(ierr);
3832f75a4aaSAlp Dener   }
3842f75a4aaSAlp Dener   PetscFunctionReturn(0);
3852f75a4aaSAlp Dener }
386c4b75bccSAlp Dener 
387c4b75bccSAlp Dener /*@C
388b6a6cedcSAlp Dener   TaoBoundSolution - Ensures that the solution vector is snapped into the bounds within a given tolerance.
389c4b75bccSAlp Dener 
390*f0fc11ceSJed Brown   Collective on Vec
391*f0fc11ceSJed Brown 
392c4b75bccSAlp Dener   Input Parameters:
393b6a6cedcSAlp Dener + X - solution vector
394b6a6cedcSAlp Dener . XL - lower bound vector
395c4b75bccSAlp Dener . XU - upper bound vector
396*f0fc11ceSJed Brown - bound_tol - absolute tolerance in enforcing the bound
397c4b75bccSAlp Dener 
398c4b75bccSAlp Dener   Output Parameters:
399*f0fc11ceSJed Brown + nDiff - total number of vector entries that have been bounded
400*f0fc11ceSJed Brown - Xout - modified solution vector satisfying bounds to bound_tol
401b6a6cedcSAlp Dener 
402b6a6cedcSAlp Dener   Level: developer
403*f0fc11ceSJed Brown 
404*f0fc11ceSJed Brown .seealso: TAOBNCG, TAOBNTL, TAOBNTR
405c4b75bccSAlp Dener @*/
4063b063059SAlp Dener PetscErrorCode TaoBoundSolution(Vec X, Vec XL, Vec XU, PetscReal bound_tol, PetscInt *nDiff, Vec Xout)
407c4b75bccSAlp Dener {
408c4b75bccSAlp Dener   PetscErrorCode    ierr;
409c4b75bccSAlp Dener   PetscInt          i,n,low,high,nDiff_loc=0;
4103b063059SAlp Dener   PetscScalar       *xout;
4113b063059SAlp Dener   const PetscScalar *x,*xl,*xu;
412c4b75bccSAlp Dener 
413c4b75bccSAlp Dener   PetscFunctionBegin;
414c4b75bccSAlp Dener   PetscValidHeaderSpecific(X,VEC_CLASSID,1);
415c4b75bccSAlp Dener   PetscValidHeaderSpecific(XL,VEC_CLASSID,2);
416c4b75bccSAlp Dener   PetscValidHeaderSpecific(XU,VEC_CLASSID,3);
4173b063059SAlp Dener   PetscValidHeaderSpecific(Xout,VEC_CLASSID,4);
418c4b75bccSAlp Dener 
419c4b75bccSAlp Dener   PetscValidType(X,1);
420c4b75bccSAlp Dener   PetscValidType(XL,2);
421c4b75bccSAlp Dener   PetscValidType(XU,3);
4223b063059SAlp Dener   PetscValidType(Xout,4);
423c4b75bccSAlp Dener   PetscCheckSameType(X,1,XL,2);
424c4b75bccSAlp Dener   PetscCheckSameType(X,1,XU,3);
4253b063059SAlp Dener   PetscCheckSameType(X,1,Xout,4);
426c4b75bccSAlp Dener   PetscCheckSameComm(X,1,XL,2);
427c4b75bccSAlp Dener   PetscCheckSameComm(X,1,XU,3);
4283b063059SAlp Dener   PetscCheckSameComm(X,1,Xout,4);
429c4b75bccSAlp Dener   VecCheckSameSize(X,1,XL,2);
430c4b75bccSAlp Dener   VecCheckSameSize(X,1,XU,3);
4313b063059SAlp Dener   VecCheckSameSize(X,1,Xout,4);
432c4b75bccSAlp Dener 
433c4b75bccSAlp Dener   ierr = VecGetOwnershipRange(X,&low,&high);CHKERRQ(ierr);
434c4b75bccSAlp Dener   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
435c4b75bccSAlp Dener   if (n>0){
43633c78596SAlp Dener     ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
43733c78596SAlp Dener     ierr = VecGetArrayRead(XL, &xl);CHKERRQ(ierr);
43833c78596SAlp Dener     ierr = VecGetArrayRead(XU, &xu);CHKERRQ(ierr);
43933c78596SAlp Dener     ierr = VecGetArray(Xout, &xout);CHKERRQ(ierr);
440c4b75bccSAlp Dener 
441c4b75bccSAlp Dener     for (i=0;i<n;++i){
44289da521bSAlp Dener       if ((xl[i] > PETSC_NINFINITY) && (x[i] <= xl[i] + bound_tol)) {
4433b063059SAlp Dener         xout[i] = xl[i]; ++nDiff_loc;
44489da521bSAlp Dener       } else if ((xu[i] < PETSC_INFINITY) && (x[i] >= xu[i] - bound_tol)) {
4453b063059SAlp Dener         xout[i] = xu[i]; ++nDiff_loc;
446c4b75bccSAlp Dener       }
447c4b75bccSAlp Dener     }
448c4b75bccSAlp Dener 
44933c78596SAlp Dener     ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);
45033c78596SAlp Dener     ierr = VecRestoreArrayRead(XL, &xl);CHKERRQ(ierr);
45133c78596SAlp Dener     ierr = VecRestoreArrayRead(XU, &xu);CHKERRQ(ierr);
45233c78596SAlp Dener     ierr = VecRestoreArray(Xout, &xout);CHKERRQ(ierr);
453c4b75bccSAlp Dener   }
4546519dc0cSAlp Dener   ierr = MPIU_Allreduce(&nDiff_loc, nDiff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)X));CHKERRQ(ierr);
455c4b75bccSAlp Dener   PetscFunctionReturn(0);
456c4b75bccSAlp Dener }
457