xref: /petsc/src/tao/bound/utils/isutil.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
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));
52*1baa6e33SBarry Smith       if (*vreduced) PetscCall(VecDestroy(vreduced));
539566063dSJacob Faibussowitsch       PetscCall(VecCreate(comm,vreduced));
549566063dSJacob Faibussowitsch       PetscCall(VecSetType(*vreduced,vtype));
55a7e14dcfSSatish Balay 
569566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(*vreduced,nreduced_local,nreduced));
579566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(*vreduced,&rlow,&rhigh));
589566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(comm,nreduced_local,rlow,1,&ident));
599566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(vfull,is,*vreduced,ident,&scatter));
609566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD));
619566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD));
629566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&scatter));
639566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&ident));
64a7e14dcfSSatish Balay       break;
65a7e14dcfSSatish Balay 
66a7e14dcfSSatish Balay     case TAO_SUBSET_MASK:
67a7e14dcfSSatish Balay     case TAO_SUBSET_MATRIXFREE:
68a7e14dcfSSatish Balay       /* vr[i] = vf[i]   if i in is
69a7e14dcfSSatish Balay        vr[i] = 0       otherwise */
706c4ed002SBarry Smith       if (!*vreduced) {
719566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(vfull,vreduced));
72a7e14dcfSSatish Balay       }
73a7e14dcfSSatish Balay 
749566063dSJacob Faibussowitsch       PetscCall(VecSet(*vreduced,maskvalue));
759566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(is,&nlocal));
769566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(vfull,&flow,&fhigh));
779566063dSJacob Faibussowitsch       PetscCall(VecGetArray(vfull,&fv));
789566063dSJacob Faibussowitsch       PetscCall(VecGetArray(*vreduced,&rv));
799566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(is,&s));
8063a3b9bcSJacob 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);
81c4b75bccSAlp Dener       for (i=0;i<nlocal;++i) {
82a7e14dcfSSatish Balay         rv[s[i]-flow] = fv[s[i]-flow];
83a7e14dcfSSatish Balay       }
849566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(is,&s));
859566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(vfull,&fv));
869566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(*vreduced,&rv));
87a7e14dcfSSatish Balay       break;
88a7e14dcfSSatish Balay     }
89a7e14dcfSSatish Balay   }
90a7e14dcfSSatish Balay   PetscFunctionReturn(0);
91a7e14dcfSSatish Balay }
92a7e14dcfSSatish Balay 
93b98f30f2SJason Sarich /*@C
94b98f30f2SJason Sarich   TaoMatGetSubMat - Gets a submatrix using the IS
95a7e14dcfSSatish Balay 
96a7e14dcfSSatish Balay   Input Parameters:
97a7e14dcfSSatish Balay + M - the full matrix (n x n)
98a7e14dcfSSatish Balay . is - the index set for the submatrix (both row and column index sets need to be the same)
99a7e14dcfSSatish Balay . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option
100b6a6cedcSAlp Dener - subset_type <TAO_SUBSET_SUBVEC,TAO_SUBSET_MASK,TAO_SUBSET_MATRIXFREE> - the method TAO is using for subsetting
101a7e14dcfSSatish Balay 
10297bb3fdcSJose E. Roman   Output Parameter:
103a7e14dcfSSatish Balay . Msub - the submatrix
104b6a6cedcSAlp Dener 
105b6a6cedcSAlp Dener   Level: developer
106a7e14dcfSSatish Balay @*/
107b98f30f2SJason Sarich PetscErrorCode TaoMatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub)
108a7e14dcfSSatish Balay {
109a7e14dcfSSatish Balay   IS             iscomp;
11062675beeSAlp Dener   PetscBool      flg = PETSC_TRUE;
11153506e15SBarry Smith 
112a7e14dcfSSatish Balay   PetscFunctionBegin;
113a7e14dcfSSatish Balay   PetscValidHeaderSpecific(M,MAT_CLASSID,1);
114a7e14dcfSSatish Balay   PetscValidHeaderSpecific(is,IS_CLASSID,2);
1159566063dSJacob Faibussowitsch   PetscCall(MatDestroy(Msub));
116a7e14dcfSSatish Balay   switch (subset_type) {
117a7e14dcfSSatish Balay   case TAO_SUBSET_SUBVEC:
1189566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub));
119a7e14dcfSSatish Balay     break;
120a7e14dcfSSatish Balay 
121a7e14dcfSSatish Balay   case TAO_SUBSET_MASK:
122a7e14dcfSSatish Balay     /* Get Reduced Hessian
123a7e14dcfSSatish Balay      Msub[i,j] = M[i,j] if i,j in Free_Local or i==j
124a7e14dcfSSatish Balay      Msub[i,j] = 0      if i!=j and i or j not in Free_Local
125a7e14dcfSSatish Balay      */
126d0609cedSBarry Smith     PetscObjectOptionsBegin((PetscObject)M);
1279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-overwrite_hessian","modify the existing hessian matrix when computing submatrices","TaoSubsetType",flg,&flg,NULL));
128d0609cedSBarry Smith     PetscOptionsEnd();
1298afaa268SBarry Smith     if (flg) {
1309566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(M, MAT_COPY_VALUES, Msub));
131a7e14dcfSSatish Balay     } else {
132a7e14dcfSSatish Balay       /* Act on hessian directly (default) */
1339566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)M));
134a7e14dcfSSatish Balay       *Msub = M;
135a7e14dcfSSatish Balay     }
136a7e14dcfSSatish Balay     /* Save the diagonal to temporary vector */
1379566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(*Msub,v1));
138a7e14dcfSSatish Balay 
139a7e14dcfSSatish Balay     /* Zero out rows and columns */
1409566063dSJacob Faibussowitsch     PetscCall(ISComplementVec(is,v1,&iscomp));
141a7e14dcfSSatish Balay 
142a7e14dcfSSatish Balay     /* Use v1 instead of 0 here because of PETSc bug */
1439566063dSJacob Faibussowitsch     PetscCall(MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1));
144a7e14dcfSSatish Balay 
1459566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iscomp));
146a7e14dcfSSatish Balay     break;
147a7e14dcfSSatish Balay   case TAO_SUBSET_MATRIXFREE:
1489566063dSJacob Faibussowitsch     PetscCall(ISComplementVec(is,v1,&iscomp));
1499566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrixFree(M,iscomp,iscomp,Msub));
1509566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iscomp));
151a7e14dcfSSatish Balay     break;
152a7e14dcfSSatish Balay   }
153a7e14dcfSSatish Balay   PetscFunctionReturn(0);
154a7e14dcfSSatish Balay }
1552f75a4aaSAlp Dener 
1562f75a4aaSAlp Dener /*@C
1572f75a4aaSAlp Dener   TaoEstimateActiveBounds - Generates index sets for variables at the lower and upper
1582f75a4aaSAlp Dener   bounds, as well as fixed variables where lower and upper bounds equal each other.
1592f75a4aaSAlp Dener 
1602f75a4aaSAlp Dener   Input Parameters:
1612f75a4aaSAlp Dener + X - solution vector
1622f75a4aaSAlp Dener . XL - lower bound vector
1632f75a4aaSAlp Dener . XU - upper bound vector
1642f75a4aaSAlp Dener . G - unprojected gradient
1650a4511e9SAlp Dener . S - step direction with which the active bounds will be estimated
166b6a6cedcSAlp Dener . W - work vector of type and size of X
1670a4511e9SAlp Dener - steplen - the step length at which the active bounds will be estimated (needs to be conservative)
1682f75a4aaSAlp Dener 
1692f75a4aaSAlp Dener   Output Parameters:
17076be6f4fSStefano Zampini + bound_tol - tolerance for the bound estimation
1712f75a4aaSAlp Dener . active_lower - index set for active variables at the lower bound
1722f75a4aaSAlp Dener . active_upper - index set for active variables at the upper bound
1732f75a4aaSAlp Dener . active_fixed - index set for fixed variables
1742f75a4aaSAlp Dener . active - index set for all active variables
175b6a6cedcSAlp Dener - inactive - complementary index set for inactive variables
1760a4511e9SAlp Dener 
1770a4511e9SAlp Dener   Notes:
1780a4511e9SAlp Dener   This estimation is based on Bertsekas' method, with a built in diagonal scaling value of 1.0e-3.
1790a4511e9SAlp Dener 
180b6a6cedcSAlp Dener   Level: developer
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   PetscReal                    wnorm;
18689da521bSAlp Dener   PetscReal                    zero = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
187c4b75bccSAlp Dener   PetscInt                     i, n_isl=0, n_isu=0, n_isf=0, n_isa=0, n_isi=0;
188c4b75bccSAlp Dener   PetscInt                     N_isl, N_isu, N_isf, N_isa, N_isi;
189c4b75bccSAlp Dener   PetscInt                     n, low, high, nDiff;
190c4b75bccSAlp Dener   PetscInt                     *isl=NULL, *isu=NULL, *isf=NULL, *isa=NULL, *isi=NULL;
1912f75a4aaSAlp Dener   const PetscScalar            *xl, *xu, *x, *g;
1926519dc0cSAlp Dener   MPI_Comm                     comm = PetscObjectComm((PetscObject)X);
1932f75a4aaSAlp Dener 
1942f75a4aaSAlp Dener   PetscFunctionBegin;
195c4b75bccSAlp Dener   PetscValidHeaderSpecific(X,VEC_CLASSID,1);
19676be6f4fSStefano Zampini   if (XL) PetscValidHeaderSpecific(XL,VEC_CLASSID,2);
19776be6f4fSStefano Zampini   if (XU) PetscValidHeaderSpecific(XU,VEC_CLASSID,3);
198c4b75bccSAlp Dener   PetscValidHeaderSpecific(G,VEC_CLASSID,4);
199c4b75bccSAlp Dener   PetscValidHeaderSpecific(S,VEC_CLASSID,5);
200c4b75bccSAlp Dener   PetscValidHeaderSpecific(W,VEC_CLASSID,6);
201c4b75bccSAlp Dener 
20276be6f4fSStefano Zampini   if (XL) PetscCheckSameType(X,1,XL,2);
20376be6f4fSStefano Zampini   if (XU) PetscCheckSameType(X,1,XU,3);
204c4b75bccSAlp Dener   PetscCheckSameType(X,1,G,4);
205c4b75bccSAlp Dener   PetscCheckSameType(X,1,S,5);
206c4b75bccSAlp Dener   PetscCheckSameType(X,1,W,6);
20776be6f4fSStefano Zampini   if (XL) PetscCheckSameComm(X,1,XL,2);
20876be6f4fSStefano Zampini   if (XU) PetscCheckSameComm(X,1,XU,3);
209c4b75bccSAlp Dener   PetscCheckSameComm(X,1,G,4);
210c4b75bccSAlp Dener   PetscCheckSameComm(X,1,S,5);
211c4b75bccSAlp Dener   PetscCheckSameComm(X,1,W,6);
21276be6f4fSStefano Zampini   if (XL) VecCheckSameSize(X,1,XL,2);
21376be6f4fSStefano Zampini   if (XU) VecCheckSameSize(X,1,XU,3);
214c4b75bccSAlp Dener   VecCheckSameSize(X,1,G,4);
215c4b75bccSAlp Dener   VecCheckSameSize(X,1,S,5);
216c4b75bccSAlp Dener   VecCheckSameSize(X,1,W,6);
217c4b75bccSAlp Dener 
2182f75a4aaSAlp Dener   /* Update the tolerance for bound detection (this is based on Bertsekas' method) */
2199566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, W));
2209566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(W, steplen, 1.0, S));
2219566063dSJacob Faibussowitsch   PetscCall(TaoBoundSolution(W, XL, XU, 0.0, &nDiff, W));
2229566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(W, 1.0, -1.0, X));
2239566063dSJacob Faibussowitsch   PetscCall(VecNorm(W, NORM_2, &wnorm));
2242f75a4aaSAlp Dener   *bound_tol = PetscMin(*bound_tol, wnorm);
2252f75a4aaSAlp Dener 
22676be6f4fSStefano Zampini   /* Clear all index sets */
22776be6f4fSStefano Zampini   PetscCall(ISDestroy(active_lower));
22876be6f4fSStefano Zampini   PetscCall(ISDestroy(active_upper));
22976be6f4fSStefano Zampini   PetscCall(ISDestroy(active_fixed));
23076be6f4fSStefano Zampini   PetscCall(ISDestroy(active));
23176be6f4fSStefano Zampini   PetscCall(ISDestroy(inactive));
23276be6f4fSStefano Zampini 
2339566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &low, &high));
2349566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
23576be6f4fSStefano Zampini   if (!XL && !XU) {
23676be6f4fSStefano Zampini     PetscCall(ISCreateStride(comm,n,low,1,inactive));
23776be6f4fSStefano Zampini     PetscFunctionReturn(0);
23876be6f4fSStefano Zampini   }
2392f75a4aaSAlp Dener   if (n>0) {
2409566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, &x));
2419566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XL, &xl));
2429566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XU, &xu));
2439566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(G, &g));
2442f75a4aaSAlp Dener 
2452f75a4aaSAlp Dener     /* Loop over variables and categorize the indexes */
2469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isl));
2479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isu));
2489566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isf));
2499566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isa));
2509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isi));
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;
25676be6f4fSStefano Zampini       } 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;
26076be6f4fSStefano Zampini       } 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 
2709566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, &x));
2719566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XL, &xl));
2729566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XU, &xu));
2739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(G, &g));
2742f75a4aaSAlp Dener   }
2752f75a4aaSAlp Dener 
276c4b75bccSAlp Dener   /* Collect global sizes */
2771c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&n_isl, &N_isl, 1, MPIU_INT, MPI_SUM, comm));
2781c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&n_isu, &N_isu, 1, MPIU_INT, MPI_SUM, comm));
2791c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&n_isf, &N_isf, 1, MPIU_INT, MPI_SUM, comm));
2801c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&n_isa, &N_isa, 1, MPIU_INT, MPI_SUM, comm));
2811c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&n_isi, &N_isi, 1, MPIU_INT, MPI_SUM, comm));
282c4b75bccSAlp Dener 
283c4b75bccSAlp Dener   /* Create index set for lower bounded variables */
284c4b75bccSAlp Dener   if (N_isl > 0) {
2859566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isl, isl, PETSC_OWN_POINTER, active_lower));
2867e9273c4SAlp Dener   } else {
2879566063dSJacob Faibussowitsch     PetscCall(PetscFree(isl));
288c4b75bccSAlp Dener   }
289c4b75bccSAlp Dener   /* Create index set for upper bounded variables */
290c4b75bccSAlp Dener   if (N_isu > 0) {
2919566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isu, isu, PETSC_OWN_POINTER, active_upper));
2927e9273c4SAlp Dener   } else {
2939566063dSJacob Faibussowitsch     PetscCall(PetscFree(isu));
294c4b75bccSAlp Dener   }
295c4b75bccSAlp Dener   /* Create index set for fixed variables */
296c4b75bccSAlp Dener   if (N_isf > 0) {
2979566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isf, isf, PETSC_OWN_POINTER, active_fixed));
2987e9273c4SAlp Dener   } else {
2999566063dSJacob Faibussowitsch     PetscCall(PetscFree(isf));
300c4b75bccSAlp Dener   }
301c4b75bccSAlp Dener   /* Create index set for all actively bounded variables */
302c4b75bccSAlp Dener   if (N_isa > 0) {
3039566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isa, isa, PETSC_OWN_POINTER, active));
3047e9273c4SAlp Dener   } else {
3059566063dSJacob Faibussowitsch     PetscCall(PetscFree(isa));
306c4b75bccSAlp Dener   }
307c4b75bccSAlp Dener   /* Create index set for all inactive variables */
308c4b75bccSAlp Dener   if (N_isi > 0) {
3099566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isi, isi, PETSC_OWN_POINTER, inactive));
3107e9273c4SAlp Dener   } else {
3119566063dSJacob Faibussowitsch     PetscCall(PetscFree(isi));
312c4b75bccSAlp Dener   }
3132f75a4aaSAlp Dener   PetscFunctionReturn(0);
3142f75a4aaSAlp Dener }
3152f75a4aaSAlp Dener 
3162f75a4aaSAlp Dener /*@C
3172f75a4aaSAlp Dener   TaoBoundStep - Ensures the correct zero or adjusted step direction
3182f75a4aaSAlp Dener   values for active variables.
3192f75a4aaSAlp Dener 
3202f75a4aaSAlp Dener   Input Parameters:
3212f75a4aaSAlp Dener + X - solution vector
3222f75a4aaSAlp Dener . XL - lower bound vector
3232f75a4aaSAlp Dener . XU - upper bound vector
3242f75a4aaSAlp Dener . active_lower - index set for lower bounded active variables
3252f75a4aaSAlp Dener . active_upper - index set for lower bounded active variables
326b6a6cedcSAlp Dener . active_fixed - index set for fixed active variables
327b6a6cedcSAlp Dener - scale - amplification factor for the step that needs to be taken on actively bounded variables
3282f75a4aaSAlp Dener 
32997bb3fdcSJose E. Roman   Output Parameter:
3302f75a4aaSAlp Dener . S - step direction to be modified
331b6a6cedcSAlp Dener 
332b6a6cedcSAlp Dener   Level: developer
3332f75a4aaSAlp Dener @*/
334c4b75bccSAlp Dener PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, PetscReal scale, Vec S)
3352f75a4aaSAlp Dener {
3362f75a4aaSAlp Dener 
3372f75a4aaSAlp Dener   Vec                          step_lower, step_upper, step_fixed;
3382f75a4aaSAlp Dener   Vec                          x_lower, x_upper;
3392f75a4aaSAlp Dener   Vec                          bound_lower, bound_upper;
3402f75a4aaSAlp Dener 
3412f75a4aaSAlp Dener   PetscFunctionBegin;
3422f75a4aaSAlp Dener   /* Adjust step for variables at the estimated lower bound */
3432f75a4aaSAlp Dener   if (active_lower) {
3449566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(S, active_lower, &step_lower));
3459566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(X, active_lower, &x_lower));
3469566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(XL, active_lower, &bound_lower));
3479566063dSJacob Faibussowitsch     PetscCall(VecCopy(bound_lower, step_lower));
3489566063dSJacob Faibussowitsch     PetscCall(VecAXPY(step_lower, -1.0, x_lower));
3499566063dSJacob Faibussowitsch     PetscCall(VecScale(step_lower, scale));
3509566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(S, active_lower, &step_lower));
3519566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(X, active_lower, &x_lower));
3529566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(XL, active_lower, &bound_lower));
3532f75a4aaSAlp Dener   }
3542f75a4aaSAlp Dener 
3552f75a4aaSAlp Dener   /* Adjust step for the variables at the estimated upper bound */
3562f75a4aaSAlp Dener   if (active_upper) {
3579566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(S, active_upper, &step_upper));
3589566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(X, active_upper, &x_upper));
3599566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(XU, active_upper, &bound_upper));
3609566063dSJacob Faibussowitsch     PetscCall(VecCopy(bound_upper, step_upper));
3619566063dSJacob Faibussowitsch     PetscCall(VecAXPY(step_upper, -1.0, x_upper));
3629566063dSJacob Faibussowitsch     PetscCall(VecScale(step_upper, scale));
3639566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(S, active_upper, &step_upper));
3649566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(X, active_upper, &x_upper));
3659566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(XU, active_upper, &bound_upper));
3662f75a4aaSAlp Dener   }
3672f75a4aaSAlp Dener 
3682f75a4aaSAlp Dener   /* Zero out step for fixed variables */
3692f75a4aaSAlp Dener   if (active_fixed) {
3709566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(S, active_fixed, &step_fixed));
3719566063dSJacob Faibussowitsch     PetscCall(VecSet(step_fixed, 0.0));
3729566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(S, active_fixed, &step_fixed));
3732f75a4aaSAlp Dener   }
3742f75a4aaSAlp Dener   PetscFunctionReturn(0);
3752f75a4aaSAlp Dener }
376c4b75bccSAlp Dener 
377c4b75bccSAlp Dener /*@C
378b6a6cedcSAlp Dener   TaoBoundSolution - Ensures that the solution vector is snapped into the bounds within a given tolerance.
379c4b75bccSAlp Dener 
380f0fc11ceSJed Brown   Collective on Vec
381f0fc11ceSJed Brown 
382c4b75bccSAlp Dener   Input Parameters:
383b6a6cedcSAlp Dener + X - solution vector
384b6a6cedcSAlp Dener . XL - lower bound vector
385c4b75bccSAlp Dener . XU - upper bound vector
386f0fc11ceSJed Brown - bound_tol - absolute tolerance in enforcing the bound
387c4b75bccSAlp Dener 
388c4b75bccSAlp Dener   Output Parameters:
389f0fc11ceSJed Brown + nDiff - total number of vector entries that have been bounded
390f0fc11ceSJed Brown - Xout - modified solution vector satisfying bounds to bound_tol
391b6a6cedcSAlp Dener 
392b6a6cedcSAlp Dener   Level: developer
393f0fc11ceSJed Brown 
394db781477SPatrick Sanan .seealso: `TAOBNCG`, `TAOBNTL`, `TAOBNTR`
395c4b75bccSAlp Dener @*/
3963b063059SAlp Dener PetscErrorCode TaoBoundSolution(Vec X, Vec XL, Vec XU, PetscReal bound_tol, PetscInt *nDiff, Vec Xout)
397c4b75bccSAlp Dener {
398c4b75bccSAlp Dener   PetscInt          i,n,low,high,nDiff_loc=0;
3993b063059SAlp Dener   PetscScalar       *xout;
4003b063059SAlp Dener   const PetscScalar *x,*xl,*xu;
401c4b75bccSAlp Dener 
402c4b75bccSAlp Dener   PetscFunctionBegin;
403c4b75bccSAlp Dener   PetscValidHeaderSpecific(X,VEC_CLASSID,1);
40476be6f4fSStefano Zampini   if (XL) PetscValidHeaderSpecific(XL,VEC_CLASSID,2);
40576be6f4fSStefano Zampini   if (XU) PetscValidHeaderSpecific(XU,VEC_CLASSID,3);
406064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Xout,VEC_CLASSID,6);
40776be6f4fSStefano Zampini   if (!XL && !XU) {
40876be6f4fSStefano Zampini     PetscCall(VecCopy(X,Xout));
40976be6f4fSStefano Zampini     *nDiff = 0.0;
41076be6f4fSStefano Zampini     PetscFunctionReturn(0);
41176be6f4fSStefano Zampini   }
412c4b75bccSAlp Dener   PetscCheckSameType(X,1,XL,2);
413c4b75bccSAlp Dener   PetscCheckSameType(X,1,XU,3);
414064a246eSJacob Faibussowitsch   PetscCheckSameType(X,1,Xout,6);
415c4b75bccSAlp Dener   PetscCheckSameComm(X,1,XL,2);
416c4b75bccSAlp Dener   PetscCheckSameComm(X,1,XU,3);
417064a246eSJacob Faibussowitsch   PetscCheckSameComm(X,1,Xout,6);
418c4b75bccSAlp Dener   VecCheckSameSize(X,1,XL,2);
419c4b75bccSAlp Dener   VecCheckSameSize(X,1,XU,3);
4203b063059SAlp Dener   VecCheckSameSize(X,1,Xout,4);
421c4b75bccSAlp Dener 
4229566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X,&low,&high));
4239566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&n));
424c4b75bccSAlp Dener   if (n>0) {
4259566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, &x));
4269566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XL, &xl));
4279566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XU, &xu));
4289566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Xout, &xout));
429c4b75bccSAlp Dener 
430c4b75bccSAlp Dener     for (i=0;i<n;++i) {
43176be6f4fSStefano Zampini       if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + bound_tol) {
4323b063059SAlp Dener         xout[i] = xl[i]; ++nDiff_loc;
43376be6f4fSStefano Zampini       } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - bound_tol) {
4343b063059SAlp Dener         xout[i] = xu[i]; ++nDiff_loc;
435c4b75bccSAlp Dener       }
436c4b75bccSAlp Dener     }
437c4b75bccSAlp Dener 
4389566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, &x));
4399566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XL, &xl));
4409566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XU, &xu));
4419566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Xout, &xout));
442c4b75bccSAlp Dener   }
4431c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&nDiff_loc, nDiff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)X)));
444c4b75bccSAlp Dener   PetscFunctionReturn(0);
445c4b75bccSAlp Dener }
446