xref: /petsc/src/tao/bound/utils/isutil.c (revision d8e47b638cf8f604a99e9678e1df24f82d959cd7)
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 
6cc4c1da9SBarry Smith /*@
73b242c63SJacob Faibussowitsch   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
123b242c63SJacob Faibussowitsch . reduced_type - the method `Tao` is using for subsetting
133b242c63SJacob Faibussowitsch - 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 
18b6a6cedcSAlp Dener   Level: developer
193b242c63SJacob Faibussowitsch 
203b242c63SJacob Faibussowitsch   Notes:
213b242c63SJacob Faibussowitsch   `maskvalue` should usually be `0.0`, unless a pointwise divide will be used.
223b242c63SJacob Faibussowitsch 
233b242c63SJacob Faibussowitsch .seealso: `TaoMatGetSubMat()`, `TaoSubsetType`
24a7e14dcfSSatish Balay @*/
TaoVecGetSubVec(Vec vfull,IS is,TaoSubsetType reduced_type,PetscReal maskvalue,Vec * vreduced)25d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoVecGetSubVec(Vec vfull, IS is, TaoSubsetType reduced_type, PetscReal maskvalue, Vec *vreduced)
26d71ae5a4SJacob Faibussowitsch {
27a7e14dcfSSatish Balay   PetscInt        nfull, nreduced, nreduced_local, rlow, rhigh, flow, fhigh;
28a7e14dcfSSatish Balay   PetscInt        i, nlocal;
29a7e14dcfSSatish Balay   PetscReal      *fv, *rv;
30a7e14dcfSSatish Balay   const PetscInt *s;
31a7e14dcfSSatish Balay   IS              ident;
32a7e14dcfSSatish Balay   VecType         vtype;
33a7e14dcfSSatish Balay   VecScatter      scatter;
34a7e14dcfSSatish Balay   MPI_Comm        comm;
35a7e14dcfSSatish Balay 
36a7e14dcfSSatish Balay   PetscFunctionBegin;
37a7e14dcfSSatish Balay   PetscValidHeaderSpecific(vfull, VEC_CLASSID, 1);
38a7e14dcfSSatish Balay   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
39a7e14dcfSSatish Balay 
409566063dSJacob Faibussowitsch   PetscCall(VecGetSize(vfull, &nfull));
419566063dSJacob Faibussowitsch   PetscCall(ISGetSize(is, &nreduced));
42a7e14dcfSSatish Balay 
43a7e14dcfSSatish Balay   if (nreduced == nfull) {
449566063dSJacob Faibussowitsch     PetscCall(VecDestroy(vreduced));
459566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(vfull, vreduced));
469566063dSJacob Faibussowitsch     PetscCall(VecCopy(vfull, *vreduced));
47a7e14dcfSSatish Balay   } else {
48a7e14dcfSSatish Balay     switch (reduced_type) {
49a7e14dcfSSatish Balay     case TAO_SUBSET_SUBVEC:
509566063dSJacob Faibussowitsch       PetscCall(VecGetType(vfull, &vtype));
519566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(vfull, &flow, &fhigh));
529566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(is, &nreduced_local));
539566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetComm((PetscObject)vfull, &comm));
541baa6e33SBarry Smith       if (*vreduced) PetscCall(VecDestroy(vreduced));
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 */
7248a46eb9SPierre Jolivet       if (!*vreduced) PetscCall(VecDuplicate(vfull, vreduced));
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);
81ad540459SPierre Jolivet       for (i = 0; i < nlocal; ++i) rv[s[i] - flow] = fv[s[i] - flow];
829566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(is, &s));
839566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(vfull, &fv));
849566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(*vreduced, &rv));
85a7e14dcfSSatish Balay       break;
86a7e14dcfSSatish Balay     }
87a7e14dcfSSatish Balay   }
883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
89a7e14dcfSSatish Balay }
90a7e14dcfSSatish Balay 
91cc4c1da9SBarry Smith /*@
922fe279fdSBarry Smith   TaoMatGetSubMat - Gets a submatrix using the `IS`
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay   Input Parameters:
953b242c63SJacob Faibussowitsch + M           - the full matrix (`n x n`)
96a7e14dcfSSatish Balay . is          - the index set for the submatrix (both row and column index sets need to be the same)
973b242c63SJacob Faibussowitsch . v1          - work vector of dimension n, needed for `TAO_SUBSET_MASK` option
983b242c63SJacob Faibussowitsch - subset_type - the method `Tao` is using for subsetting
99a7e14dcfSSatish Balay 
10097bb3fdcSJose E. Roman   Output Parameter:
101a7e14dcfSSatish Balay . Msub - the submatrix
102b6a6cedcSAlp Dener 
103b6a6cedcSAlp Dener   Level: developer
1043b242c63SJacob Faibussowitsch 
1053b242c63SJacob Faibussowitsch .seealso: `TaoVecGetSubVec()`, `TaoSubsetType`
106a7e14dcfSSatish Balay @*/
TaoMatGetSubMat(Mat M,IS is,Vec v1,TaoSubsetType subset_type,Mat * Msub)107d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoMatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub)
108d71ae5a4SJacob Faibussowitsch {
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) {
117d71ae5a4SJacob Faibussowitsch   case TAO_SUBSET_SUBVEC:
118d71ae5a4SJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub));
119d71ae5a4SJacob Faibussowitsch     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   }
1533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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
1663b242c63SJacob Faibussowitsch . 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 
177b6a6cedcSAlp Dener   Level: developer
1783b242c63SJacob Faibussowitsch 
1793b242c63SJacob Faibussowitsch   Notes:
1803b242c63SJacob Faibussowitsch   This estimation is based on Bertsekas' method, with a built in diagonal scaling value of `1.0e-3`.
1813b242c63SJacob Faibussowitsch 
1823b242c63SJacob Faibussowitsch .seealso: `TAOBNCG`, `TAOBNTL`, `TAOBNTR`, `TaoBoundSolution()`
1832f75a4aaSAlp Dener @*/
TaoEstimateActiveBounds(Vec X,Vec XL,Vec XU,Vec G,Vec S,Vec W,PetscReal steplen,PetscReal * bound_tol,IS * active_lower,IS * active_upper,IS * active_fixed,IS * active,IS * inactive)184d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoEstimateActiveBounds(Vec X, Vec XL, Vec XU, Vec G, Vec S, Vec W, PetscReal steplen, PetscReal *bound_tol, IS *active_lower, IS *active_upper, IS *active_fixed, IS *active, IS *inactive)
185d71ae5a4SJacob Faibussowitsch {
1862f75a4aaSAlp Dener   PetscReal          wnorm;
18789da521bSAlp Dener   PetscReal          zero = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0);
188c4b75bccSAlp Dener   PetscInt           i, n_isl = 0, n_isu = 0, n_isf = 0, n_isa = 0, n_isi = 0;
189c4b75bccSAlp Dener   PetscInt           N_isl, N_isu, N_isf, N_isa, N_isi;
190c4b75bccSAlp Dener   PetscInt           n, low, high, nDiff;
191c4b75bccSAlp Dener   PetscInt          *isl = NULL, *isu = NULL, *isf = NULL, *isa = NULL, *isi = NULL;
1922f75a4aaSAlp Dener   const PetscScalar *xl, *xu, *x, *g;
1936519dc0cSAlp Dener   MPI_Comm           comm = PetscObjectComm((PetscObject)X);
1942f75a4aaSAlp Dener 
1952f75a4aaSAlp Dener   PetscFunctionBegin;
196c4b75bccSAlp Dener   PetscValidHeaderSpecific(X, VEC_CLASSID, 1);
19776be6f4fSStefano Zampini   if (XL) PetscValidHeaderSpecific(XL, VEC_CLASSID, 2);
19876be6f4fSStefano Zampini   if (XU) PetscValidHeaderSpecific(XU, VEC_CLASSID, 3);
199c4b75bccSAlp Dener   PetscValidHeaderSpecific(G, VEC_CLASSID, 4);
200c4b75bccSAlp Dener   PetscValidHeaderSpecific(S, VEC_CLASSID, 5);
201c4b75bccSAlp Dener   PetscValidHeaderSpecific(W, VEC_CLASSID, 6);
202c4b75bccSAlp Dener 
20376be6f4fSStefano Zampini   if (XL) PetscCheckSameType(X, 1, XL, 2);
20476be6f4fSStefano Zampini   if (XU) PetscCheckSameType(X, 1, XU, 3);
205c4b75bccSAlp Dener   PetscCheckSameType(X, 1, G, 4);
206c4b75bccSAlp Dener   PetscCheckSameType(X, 1, S, 5);
207c4b75bccSAlp Dener   PetscCheckSameType(X, 1, W, 6);
20876be6f4fSStefano Zampini   if (XL) PetscCheckSameComm(X, 1, XL, 2);
20976be6f4fSStefano Zampini   if (XU) PetscCheckSameComm(X, 1, XU, 3);
210c4b75bccSAlp Dener   PetscCheckSameComm(X, 1, G, 4);
211c4b75bccSAlp Dener   PetscCheckSameComm(X, 1, S, 5);
212c4b75bccSAlp Dener   PetscCheckSameComm(X, 1, W, 6);
21376be6f4fSStefano Zampini   if (XL) VecCheckSameSize(X, 1, XL, 2);
21476be6f4fSStefano Zampini   if (XU) VecCheckSameSize(X, 1, XU, 3);
215c4b75bccSAlp Dener   VecCheckSameSize(X, 1, G, 4);
216c4b75bccSAlp Dener   VecCheckSameSize(X, 1, S, 5);
217c4b75bccSAlp Dener   VecCheckSameSize(X, 1, W, 6);
218c4b75bccSAlp Dener 
2192f75a4aaSAlp Dener   /* Update the tolerance for bound detection (this is based on Bertsekas' method) */
2209566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, W));
2219566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(W, steplen, 1.0, S));
2229566063dSJacob Faibussowitsch   PetscCall(TaoBoundSolution(W, XL, XU, 0.0, &nDiff, W));
2239566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(W, 1.0, -1.0, X));
2249566063dSJacob Faibussowitsch   PetscCall(VecNorm(W, NORM_2, &wnorm));
2252f75a4aaSAlp Dener   *bound_tol = PetscMin(*bound_tol, wnorm);
2262f75a4aaSAlp Dener 
22776be6f4fSStefano Zampini   /* Clear all index sets */
22876be6f4fSStefano Zampini   PetscCall(ISDestroy(active_lower));
22976be6f4fSStefano Zampini   PetscCall(ISDestroy(active_upper));
23076be6f4fSStefano Zampini   PetscCall(ISDestroy(active_fixed));
23176be6f4fSStefano Zampini   PetscCall(ISDestroy(active));
23276be6f4fSStefano Zampini   PetscCall(ISDestroy(inactive));
23376be6f4fSStefano Zampini 
2349566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &low, &high));
2359566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
23676be6f4fSStefano Zampini   if (!XL && !XU) {
23776be6f4fSStefano Zampini     PetscCall(ISCreateStride(comm, n, low, 1, inactive));
2383ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
23976be6f4fSStefano Zampini   }
2402f75a4aaSAlp Dener   if (n > 0) {
2419566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, &x));
2429566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XL, &xl));
2439566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XU, &xu));
2449566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(G, &g));
2452f75a4aaSAlp Dener 
2462f75a4aaSAlp Dener     /* Loop over variables and categorize the indexes */
2479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isl));
2489566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isu));
2499566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isf));
2509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isa));
2519566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &isi));
252c4b75bccSAlp Dener     for (i = 0; i < n; ++i) {
25366b4d56fSTodd Munson       if (xl[i] == xu[i]) {
254c4b75bccSAlp Dener         /* Fixed variables */
2559371c9d4SSatish Balay         isf[n_isf] = low + i;
2569371c9d4SSatish Balay         ++n_isf;
2579371c9d4SSatish Balay         isa[n_isa] = low + i;
2589371c9d4SSatish Balay         ++n_isa;
25976be6f4fSStefano Zampini       } else if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + *bound_tol && g[i] > zero) {
260c4b75bccSAlp Dener         /* Lower bounded variables */
2619371c9d4SSatish Balay         isl[n_isl] = low + i;
2629371c9d4SSatish Balay         ++n_isl;
2639371c9d4SSatish Balay         isa[n_isa] = low + i;
2649371c9d4SSatish Balay         ++n_isa;
26576be6f4fSStefano Zampini       } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - *bound_tol && g[i] < zero) {
266c4b75bccSAlp Dener         /* Upper bounded variables */
2679371c9d4SSatish Balay         isu[n_isu] = low + i;
2689371c9d4SSatish Balay         ++n_isu;
2699371c9d4SSatish Balay         isa[n_isa] = low + i;
2709371c9d4SSatish Balay         ++n_isa;
271c4b75bccSAlp Dener       } else {
272c4b75bccSAlp Dener         /* Inactive variables */
2739371c9d4SSatish Balay         isi[n_isi] = low + i;
2749371c9d4SSatish Balay         ++n_isi;
2752f75a4aaSAlp Dener       }
2762f75a4aaSAlp Dener     }
2772f75a4aaSAlp Dener 
2789566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, &x));
2799566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XL, &xl));
2809566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XU, &xu));
2819566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(G, &g));
2822f75a4aaSAlp Dener   }
2832f75a4aaSAlp Dener 
284c4b75bccSAlp Dener   /* Collect global sizes */
285*462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&n_isl, &N_isl, 1, MPIU_INT, MPI_SUM, comm));
286*462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&n_isu, &N_isu, 1, MPIU_INT, MPI_SUM, comm));
287*462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&n_isf, &N_isf, 1, MPIU_INT, MPI_SUM, comm));
288*462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&n_isa, &N_isa, 1, MPIU_INT, MPI_SUM, comm));
289*462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&n_isi, &N_isi, 1, MPIU_INT, MPI_SUM, comm));
290c4b75bccSAlp Dener 
291c4b75bccSAlp Dener   /* Create index set for lower bounded variables */
292c4b75bccSAlp Dener   if (N_isl > 0) {
2939566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isl, isl, PETSC_OWN_POINTER, active_lower));
2947e9273c4SAlp Dener   } else {
2959566063dSJacob Faibussowitsch     PetscCall(PetscFree(isl));
296c4b75bccSAlp Dener   }
297c4b75bccSAlp Dener   /* Create index set for upper bounded variables */
298c4b75bccSAlp Dener   if (N_isu > 0) {
2999566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isu, isu, PETSC_OWN_POINTER, active_upper));
3007e9273c4SAlp Dener   } else {
3019566063dSJacob Faibussowitsch     PetscCall(PetscFree(isu));
302c4b75bccSAlp Dener   }
303c4b75bccSAlp Dener   /* Create index set for fixed variables */
304c4b75bccSAlp Dener   if (N_isf > 0) {
3059566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isf, isf, PETSC_OWN_POINTER, active_fixed));
3067e9273c4SAlp Dener   } else {
3079566063dSJacob Faibussowitsch     PetscCall(PetscFree(isf));
308c4b75bccSAlp Dener   }
309c4b75bccSAlp Dener   /* Create index set for all actively bounded variables */
310c4b75bccSAlp Dener   if (N_isa > 0) {
3119566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isa, isa, PETSC_OWN_POINTER, active));
3127e9273c4SAlp Dener   } else {
3139566063dSJacob Faibussowitsch     PetscCall(PetscFree(isa));
314c4b75bccSAlp Dener   }
315c4b75bccSAlp Dener   /* Create index set for all inactive variables */
316c4b75bccSAlp Dener   if (N_isi > 0) {
3179566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_isi, isi, PETSC_OWN_POINTER, inactive));
3187e9273c4SAlp Dener   } else {
3199566063dSJacob Faibussowitsch     PetscCall(PetscFree(isi));
320c4b75bccSAlp Dener   }
3213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3222f75a4aaSAlp Dener }
3232f75a4aaSAlp Dener 
324cc4c1da9SBarry Smith /*@
3253b242c63SJacob Faibussowitsch   TaoBoundStep - Ensures the correct zero or adjusted step direction values for active
3263b242c63SJacob Faibussowitsch   variables.
3272f75a4aaSAlp Dener 
3282f75a4aaSAlp Dener   Input Parameters:
3292f75a4aaSAlp Dener + X            - solution vector
3302f75a4aaSAlp Dener . XL           - lower bound vector
3312f75a4aaSAlp Dener . XU           - upper bound vector
3322f75a4aaSAlp Dener . active_lower - index set for lower bounded active variables
3332f75a4aaSAlp Dener . active_upper - index set for lower bounded active variables
334b6a6cedcSAlp Dener . active_fixed - index set for fixed active variables
335b6a6cedcSAlp Dener - scale        - amplification factor for the step that needs to be taken on actively bounded variables
3362f75a4aaSAlp Dener 
33797bb3fdcSJose E. Roman   Output Parameter:
3382f75a4aaSAlp Dener . S - step direction to be modified
339b6a6cedcSAlp Dener 
340b6a6cedcSAlp Dener   Level: developer
3413b242c63SJacob Faibussowitsch 
3423b242c63SJacob Faibussowitsch .seealso: `TAOBNCG`, `TAOBNTL`, `TAOBNTR`, `TaoBoundSolution()`
3432f75a4aaSAlp Dener @*/
TaoBoundStep(Vec X,Vec XL,Vec XU,IS active_lower,IS active_upper,IS active_fixed,PetscReal scale,Vec S)344d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, PetscReal scale, Vec S)
345d71ae5a4SJacob Faibussowitsch {
3462f75a4aaSAlp Dener   Vec step_lower, step_upper, step_fixed;
3472f75a4aaSAlp Dener   Vec x_lower, x_upper;
3482f75a4aaSAlp Dener   Vec bound_lower, bound_upper;
3492f75a4aaSAlp Dener 
3502f75a4aaSAlp Dener   PetscFunctionBegin;
3512f75a4aaSAlp Dener   /* Adjust step for variables at the estimated lower bound */
3522f75a4aaSAlp Dener   if (active_lower) {
3539566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(S, active_lower, &step_lower));
3549566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(X, active_lower, &x_lower));
3559566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(XL, active_lower, &bound_lower));
3569566063dSJacob Faibussowitsch     PetscCall(VecCopy(bound_lower, step_lower));
3579566063dSJacob Faibussowitsch     PetscCall(VecAXPY(step_lower, -1.0, x_lower));
3589566063dSJacob Faibussowitsch     PetscCall(VecScale(step_lower, scale));
3599566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(S, active_lower, &step_lower));
3609566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(X, active_lower, &x_lower));
3619566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(XL, active_lower, &bound_lower));
3622f75a4aaSAlp Dener   }
3632f75a4aaSAlp Dener 
3642f75a4aaSAlp Dener   /* Adjust step for the variables at the estimated upper bound */
3652f75a4aaSAlp Dener   if (active_upper) {
3669566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(S, active_upper, &step_upper));
3679566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(X, active_upper, &x_upper));
3689566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(XU, active_upper, &bound_upper));
3699566063dSJacob Faibussowitsch     PetscCall(VecCopy(bound_upper, step_upper));
3709566063dSJacob Faibussowitsch     PetscCall(VecAXPY(step_upper, -1.0, x_upper));
3719566063dSJacob Faibussowitsch     PetscCall(VecScale(step_upper, scale));
3729566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(S, active_upper, &step_upper));
3739566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(X, active_upper, &x_upper));
3749566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(XU, active_upper, &bound_upper));
3752f75a4aaSAlp Dener   }
3762f75a4aaSAlp Dener 
3772f75a4aaSAlp Dener   /* Zero out step for fixed variables */
3782f75a4aaSAlp Dener   if (active_fixed) {
3799566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(S, active_fixed, &step_fixed));
3809566063dSJacob Faibussowitsch     PetscCall(VecSet(step_fixed, 0.0));
3819566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(S, active_fixed, &step_fixed));
3822f75a4aaSAlp Dener   }
3833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3842f75a4aaSAlp Dener }
385c4b75bccSAlp Dener 
386cc4c1da9SBarry Smith /*@
387b6a6cedcSAlp Dener   TaoBoundSolution - Ensures that the solution vector is snapped into the bounds within a given tolerance.
388c4b75bccSAlp Dener 
389c3339decSBarry Smith   Collective
390f0fc11ceSJed Brown 
391c4b75bccSAlp Dener   Input Parameters:
392b6a6cedcSAlp Dener + X         - solution vector
393b6a6cedcSAlp Dener . XL        - lower bound vector
394c4b75bccSAlp Dener . XU        - upper bound vector
395f0fc11ceSJed Brown - bound_tol - absolute tolerance in enforcing the bound
396c4b75bccSAlp Dener 
397c4b75bccSAlp Dener   Output Parameters:
398f0fc11ceSJed Brown + nDiff - total number of vector entries that have been bounded
3993b242c63SJacob Faibussowitsch - Xout  - modified solution vector satisfying bounds to `bound_tol`
400b6a6cedcSAlp Dener 
401b6a6cedcSAlp Dener   Level: developer
402f0fc11ceSJed Brown 
4033b242c63SJacob Faibussowitsch .seealso: `TAOBNCG`, `TAOBNTL`, `TAOBNTR`, `TaoBoundStep()`
404c4b75bccSAlp Dener @*/
TaoBoundSolution(Vec X,Vec XL,Vec XU,PetscReal bound_tol,PetscInt * nDiff,Vec Xout)405d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBoundSolution(Vec X, Vec XL, Vec XU, PetscReal bound_tol, PetscInt *nDiff, Vec Xout)
406d71ae5a4SJacob Faibussowitsch {
407c4b75bccSAlp Dener   PetscInt           i, n, low, high, nDiff_loc = 0;
4083b063059SAlp Dener   PetscScalar       *xout;
4093b063059SAlp Dener   const PetscScalar *x, *xl, *xu;
410c4b75bccSAlp Dener 
411c4b75bccSAlp Dener   PetscFunctionBegin;
412c4b75bccSAlp Dener   PetscValidHeaderSpecific(X, VEC_CLASSID, 1);
41376be6f4fSStefano Zampini   if (XL) PetscValidHeaderSpecific(XL, VEC_CLASSID, 2);
41476be6f4fSStefano Zampini   if (XU) PetscValidHeaderSpecific(XU, VEC_CLASSID, 3);
415064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(Xout, VEC_CLASSID, 6);
41676be6f4fSStefano Zampini   if (!XL && !XU) {
41776be6f4fSStefano Zampini     PetscCall(VecCopy(X, Xout));
41876be6f4fSStefano Zampini     *nDiff = 0.0;
4193ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
42076be6f4fSStefano Zampini   }
421c4b75bccSAlp Dener   PetscCheckSameType(X, 1, XL, 2);
422c4b75bccSAlp Dener   PetscCheckSameType(X, 1, XU, 3);
423064a246eSJacob Faibussowitsch   PetscCheckSameType(X, 1, Xout, 6);
424c4b75bccSAlp Dener   PetscCheckSameComm(X, 1, XL, 2);
425c4b75bccSAlp Dener   PetscCheckSameComm(X, 1, XU, 3);
426064a246eSJacob Faibussowitsch   PetscCheckSameComm(X, 1, Xout, 6);
427c4b75bccSAlp Dener   VecCheckSameSize(X, 1, XL, 2);
428c4b75bccSAlp Dener   VecCheckSameSize(X, 1, XU, 3);
4293b063059SAlp Dener   VecCheckSameSize(X, 1, Xout, 4);
430c4b75bccSAlp Dener 
4319566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &low, &high));
4329566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
433c4b75bccSAlp Dener   if (n > 0) {
4349566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, &x));
4359566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XL, &xl));
4369566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(XU, &xu));
4379566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Xout, &xout));
438c4b75bccSAlp Dener 
439c4b75bccSAlp Dener     for (i = 0; i < n; ++i) {
44076be6f4fSStefano Zampini       if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + bound_tol) {
4419371c9d4SSatish Balay         xout[i] = xl[i];
4429371c9d4SSatish Balay         ++nDiff_loc;
44376be6f4fSStefano Zampini       } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - bound_tol) {
4449371c9d4SSatish Balay         xout[i] = xu[i];
4459371c9d4SSatish Balay         ++nDiff_loc;
446c4b75bccSAlp Dener       }
447c4b75bccSAlp Dener     }
448c4b75bccSAlp Dener 
4499566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, &x));
4509566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XL, &xl));
4519566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(XU, &xu));
4529566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Xout, &xout));
453c4b75bccSAlp Dener   }
454*462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&nDiff_loc, nDiff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)X)));
4553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
456c4b75bccSAlp Dener }
457