xref: /petsc/src/tao/bound/utils/isutil.c (revision 76be6f4ff3bd4e251c19fc00ebbebfd58b6e7589)
1 #include <petsctao.h> /*I "petsctao.h" I*/
2 #include <petsc/private/vecimpl.h>
3 #include <petsc/private/taoimpl.h>
4 #include <../src/tao/matrix/submatfree.h>
5 
6 /*@C
7   TaoVecGetSubVec - Gets a subvector using the IS
8 
9   Input Parameters:
10 + vfull - the full matrix
11 . is - the index set for the subvector
12 . reduced_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK,  TAO_SUBSET_MATRIXFREE)
13 - maskvalue - the value to set the unused vector elements to (for TAO_SUBSET_MASK or TAO_SUBSET_MATRIXFREE)
14 
15   Output Parameter:
16 . vreduced - the subvector
17 
18   Notes:
19   maskvalue should usually be 0.0, unless a pointwise divide will be used.
20 
21   Level: developer
22 @*/
23 PetscErrorCode TaoVecGetSubVec(Vec vfull, IS is, TaoSubsetType reduced_type, PetscReal maskvalue, Vec *vreduced)
24 {
25   PetscInt       nfull,nreduced,nreduced_local,rlow,rhigh,flow,fhigh;
26   PetscInt       i,nlocal;
27   PetscReal      *fv,*rv;
28   const PetscInt *s;
29   IS             ident;
30   VecType        vtype;
31   VecScatter     scatter;
32   MPI_Comm       comm;
33 
34   PetscFunctionBegin;
35   PetscValidHeaderSpecific(vfull,VEC_CLASSID,1);
36   PetscValidHeaderSpecific(is,IS_CLASSID,2);
37 
38   PetscCall(VecGetSize(vfull, &nfull));
39   PetscCall(ISGetSize(is, &nreduced));
40 
41   if (nreduced == nfull) {
42     PetscCall(VecDestroy(vreduced));
43     PetscCall(VecDuplicate(vfull,vreduced));
44     PetscCall(VecCopy(vfull,*vreduced));
45   } else {
46     switch (reduced_type) {
47     case TAO_SUBSET_SUBVEC:
48       PetscCall(VecGetType(vfull,&vtype));
49       PetscCall(VecGetOwnershipRange(vfull,&flow,&fhigh));
50       PetscCall(ISGetLocalSize(is,&nreduced_local));
51       PetscCall(PetscObjectGetComm((PetscObject)vfull,&comm));
52       if (*vreduced) {
53         PetscCall(VecDestroy(vreduced));
54       }
55       PetscCall(VecCreate(comm,vreduced));
56       PetscCall(VecSetType(*vreduced,vtype));
57 
58       PetscCall(VecSetSizes(*vreduced,nreduced_local,nreduced));
59       PetscCall(VecGetOwnershipRange(*vreduced,&rlow,&rhigh));
60       PetscCall(ISCreateStride(comm,nreduced_local,rlow,1,&ident));
61       PetscCall(VecScatterCreate(vfull,is,*vreduced,ident,&scatter));
62       PetscCall(VecScatterBegin(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD));
63       PetscCall(VecScatterEnd(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD));
64       PetscCall(VecScatterDestroy(&scatter));
65       PetscCall(ISDestroy(&ident));
66       break;
67 
68     case TAO_SUBSET_MASK:
69     case TAO_SUBSET_MATRIXFREE:
70       /* vr[i] = vf[i]   if i in is
71        vr[i] = 0       otherwise */
72       if (!*vreduced) {
73         PetscCall(VecDuplicate(vfull,vreduced));
74       }
75 
76       PetscCall(VecSet(*vreduced,maskvalue));
77       PetscCall(ISGetLocalSize(is,&nlocal));
78       PetscCall(VecGetOwnershipRange(vfull,&flow,&fhigh));
79       PetscCall(VecGetArray(vfull,&fv));
80       PetscCall(VecGetArray(*vreduced,&rv));
81       PetscCall(ISGetIndices(is,&s));
82       PetscCheck(nlocal <= (fhigh-flow),PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"IS local size %" PetscInt_FMT " > Vec local size %" PetscInt_FMT,nlocal,fhigh-flow);
83       for (i=0;i<nlocal;++i) {
84         rv[s[i]-flow] = fv[s[i]-flow];
85       }
86       PetscCall(ISRestoreIndices(is,&s));
87       PetscCall(VecRestoreArray(vfull,&fv));
88       PetscCall(VecRestoreArray(*vreduced,&rv));
89       break;
90     }
91   }
92   PetscFunctionReturn(0);
93 }
94 
95 /*@C
96   TaoMatGetSubMat - Gets a submatrix using the IS
97 
98   Input Parameters:
99 + M - the full matrix (n x n)
100 . is - the index set for the submatrix (both row and column index sets need to be the same)
101 . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option
102 - subset_type <TAO_SUBSET_SUBVEC,TAO_SUBSET_MASK,TAO_SUBSET_MATRIXFREE> - the method TAO is using for subsetting
103 
104   Output Parameter:
105 . Msub - the submatrix
106 
107   Level: developer
108 @*/
109 PetscErrorCode TaoMatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub)
110 {
111   IS             iscomp;
112   PetscBool      flg = PETSC_TRUE;
113 
114   PetscFunctionBegin;
115   PetscValidHeaderSpecific(M,MAT_CLASSID,1);
116   PetscValidHeaderSpecific(is,IS_CLASSID,2);
117   PetscCall(MatDestroy(Msub));
118   switch (subset_type) {
119   case TAO_SUBSET_SUBVEC:
120     PetscCall(MatCreateSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub));
121     break;
122 
123   case TAO_SUBSET_MASK:
124     /* Get Reduced Hessian
125      Msub[i,j] = M[i,j] if i,j in Free_Local or i==j
126      Msub[i,j] = 0      if i!=j and i or j not in Free_Local
127      */
128     PetscObjectOptionsBegin((PetscObject)M);
129     PetscCall(PetscOptionsBool("-overwrite_hessian","modify the existing hessian matrix when computing submatrices","TaoSubsetType",flg,&flg,NULL));
130     PetscOptionsEnd();
131     if (flg) {
132       PetscCall(MatDuplicate(M, MAT_COPY_VALUES, Msub));
133     } else {
134       /* Act on hessian directly (default) */
135       PetscCall(PetscObjectReference((PetscObject)M));
136       *Msub = M;
137     }
138     /* Save the diagonal to temporary vector */
139     PetscCall(MatGetDiagonal(*Msub,v1));
140 
141     /* Zero out rows and columns */
142     PetscCall(ISComplementVec(is,v1,&iscomp));
143 
144     /* Use v1 instead of 0 here because of PETSc bug */
145     PetscCall(MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1));
146 
147     PetscCall(ISDestroy(&iscomp));
148     break;
149   case TAO_SUBSET_MATRIXFREE:
150     PetscCall(ISComplementVec(is,v1,&iscomp));
151     PetscCall(MatCreateSubMatrixFree(M,iscomp,iscomp,Msub));
152     PetscCall(ISDestroy(&iscomp));
153     break;
154   }
155   PetscFunctionReturn(0);
156 }
157 
158 /*@C
159   TaoEstimateActiveBounds - Generates index sets for variables at the lower and upper
160   bounds, as well as fixed variables where lower and upper bounds equal each other.
161 
162   Input Parameters:
163 + X - solution vector
164 . XL - lower bound vector
165 . XU - upper bound vector
166 . G - unprojected gradient
167 . S - step direction with which the active bounds will be estimated
168 . W - work vector of type and size of X
169 - steplen - the step length at which the active bounds will be estimated (needs to be conservative)
170 
171   Output Parameters:
172 + bound_tol - tolerance for the bound estimation
173 . active_lower - index set for active variables at the lower bound
174 . active_upper - index set for active variables at the upper bound
175 . active_fixed - index set for fixed variables
176 . active - index set for all active variables
177 - inactive - complementary index set for inactive variables
178 
179   Notes:
180   This estimation is based on Bertsekas' method, with a built in diagonal scaling value of 1.0e-3.
181 
182   Level: developer
183 @*/
184 PetscErrorCode TaoEstimateActiveBounds(Vec X, Vec XL, Vec XU, Vec G, Vec S, Vec W, PetscReal steplen, PetscReal *bound_tol,
185                                        IS *active_lower, IS *active_upper, IS *active_fixed, IS *active, IS *inactive)
186 {
187   PetscReal                    wnorm;
188   PetscReal                    zero = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
189   PetscInt                     i, n_isl=0, n_isu=0, n_isf=0, n_isa=0, n_isi=0;
190   PetscInt                     N_isl, N_isu, N_isf, N_isa, N_isi;
191   PetscInt                     n, low, high, nDiff;
192   PetscInt                     *isl=NULL, *isu=NULL, *isf=NULL, *isa=NULL, *isi=NULL;
193   const PetscScalar            *xl, *xu, *x, *g;
194   MPI_Comm                     comm = PetscObjectComm((PetscObject)X);
195 
196   PetscFunctionBegin;
197   PetscValidHeaderSpecific(X,VEC_CLASSID,1);
198   if (XL) PetscValidHeaderSpecific(XL,VEC_CLASSID,2);
199   if (XU) PetscValidHeaderSpecific(XU,VEC_CLASSID,3);
200   PetscValidHeaderSpecific(G,VEC_CLASSID,4);
201   PetscValidHeaderSpecific(S,VEC_CLASSID,5);
202   PetscValidHeaderSpecific(W,VEC_CLASSID,6);
203 
204   if (XL) PetscCheckSameType(X,1,XL,2);
205   if (XU) PetscCheckSameType(X,1,XU,3);
206   PetscCheckSameType(X,1,G,4);
207   PetscCheckSameType(X,1,S,5);
208   PetscCheckSameType(X,1,W,6);
209   if (XL) PetscCheckSameComm(X,1,XL,2);
210   if (XU) PetscCheckSameComm(X,1,XU,3);
211   PetscCheckSameComm(X,1,G,4);
212   PetscCheckSameComm(X,1,S,5);
213   PetscCheckSameComm(X,1,W,6);
214   if (XL) VecCheckSameSize(X,1,XL,2);
215   if (XU) VecCheckSameSize(X,1,XU,3);
216   VecCheckSameSize(X,1,G,4);
217   VecCheckSameSize(X,1,S,5);
218   VecCheckSameSize(X,1,W,6);
219 
220   /* Update the tolerance for bound detection (this is based on Bertsekas' method) */
221   PetscCall(VecCopy(X, W));
222   PetscCall(VecAXPBY(W, steplen, 1.0, S));
223   PetscCall(TaoBoundSolution(W, XL, XU, 0.0, &nDiff, W));
224   PetscCall(VecAXPBY(W, 1.0, -1.0, X));
225   PetscCall(VecNorm(W, NORM_2, &wnorm));
226   *bound_tol = PetscMin(*bound_tol, wnorm);
227 
228   /* Clear all index sets */
229   PetscCall(ISDestroy(active_lower));
230   PetscCall(ISDestroy(active_upper));
231   PetscCall(ISDestroy(active_fixed));
232   PetscCall(ISDestroy(active));
233   PetscCall(ISDestroy(inactive));
234 
235   PetscCall(VecGetOwnershipRange(X, &low, &high));
236   PetscCall(VecGetLocalSize(X, &n));
237   if (!XL && !XU) {
238     PetscCall(ISCreateStride(comm,n,low,1,inactive));
239     PetscFunctionReturn(0);
240   }
241   if (n>0) {
242     PetscCall(VecGetArrayRead(X, &x));
243     PetscCall(VecGetArrayRead(XL, &xl));
244     PetscCall(VecGetArrayRead(XU, &xu));
245     PetscCall(VecGetArrayRead(G, &g));
246 
247     /* Loop over variables and categorize the indexes */
248     PetscCall(PetscMalloc1(n, &isl));
249     PetscCall(PetscMalloc1(n, &isu));
250     PetscCall(PetscMalloc1(n, &isf));
251     PetscCall(PetscMalloc1(n, &isa));
252     PetscCall(PetscMalloc1(n, &isi));
253     for (i=0; i<n; ++i) {
254       if (xl[i] == xu[i]) {
255         /* Fixed variables */
256         isf[n_isf]=low+i; ++n_isf;
257         isa[n_isa]=low+i; ++n_isa;
258       } else if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + *bound_tol && g[i] > zero) {
259         /* Lower bounded variables */
260         isl[n_isl]=low+i; ++n_isl;
261         isa[n_isa]=low+i; ++n_isa;
262       } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - *bound_tol && g[i] < zero) {
263         /* Upper bounded variables */
264         isu[n_isu]=low+i; ++n_isu;
265         isa[n_isa]=low+i; ++n_isa;
266       } else {
267         /* Inactive variables */
268         isi[n_isi]=low+i; ++n_isi;
269       }
270     }
271 
272     PetscCall(VecRestoreArrayRead(X, &x));
273     PetscCall(VecRestoreArrayRead(XL, &xl));
274     PetscCall(VecRestoreArrayRead(XU, &xu));
275     PetscCall(VecRestoreArrayRead(G, &g));
276   }
277 
278   /* Collect global sizes */
279   PetscCall(MPIU_Allreduce(&n_isl, &N_isl, 1, MPIU_INT, MPI_SUM, comm));
280   PetscCall(MPIU_Allreduce(&n_isu, &N_isu, 1, MPIU_INT, MPI_SUM, comm));
281   PetscCall(MPIU_Allreduce(&n_isf, &N_isf, 1, MPIU_INT, MPI_SUM, comm));
282   PetscCall(MPIU_Allreduce(&n_isa, &N_isa, 1, MPIU_INT, MPI_SUM, comm));
283   PetscCall(MPIU_Allreduce(&n_isi, &N_isi, 1, MPIU_INT, MPI_SUM, comm));
284 
285   /* Create index set for lower bounded variables */
286   if (N_isl > 0) {
287     PetscCall(ISCreateGeneral(comm, n_isl, isl, PETSC_OWN_POINTER, active_lower));
288   } else {
289     PetscCall(PetscFree(isl));
290   }
291   /* Create index set for upper bounded variables */
292   if (N_isu > 0) {
293     PetscCall(ISCreateGeneral(comm, n_isu, isu, PETSC_OWN_POINTER, active_upper));
294   } else {
295     PetscCall(PetscFree(isu));
296   }
297   /* Create index set for fixed variables */
298   if (N_isf > 0) {
299     PetscCall(ISCreateGeneral(comm, n_isf, isf, PETSC_OWN_POINTER, active_fixed));
300   } else {
301     PetscCall(PetscFree(isf));
302   }
303   /* Create index set for all actively bounded variables */
304   if (N_isa > 0) {
305     PetscCall(ISCreateGeneral(comm, n_isa, isa, PETSC_OWN_POINTER, active));
306   } else {
307     PetscCall(PetscFree(isa));
308   }
309   /* Create index set for all inactive variables */
310   if (N_isi > 0) {
311     PetscCall(ISCreateGeneral(comm, n_isi, isi, PETSC_OWN_POINTER, inactive));
312   } else {
313     PetscCall(PetscFree(isi));
314   }
315   PetscFunctionReturn(0);
316 }
317 
318 /*@C
319   TaoBoundStep - Ensures the correct zero or adjusted step direction
320   values for active variables.
321 
322   Input Parameters:
323 + X - solution vector
324 . XL - lower bound vector
325 . XU - upper bound vector
326 . active_lower - index set for lower bounded active variables
327 . active_upper - index set for lower bounded active variables
328 . active_fixed - index set for fixed active variables
329 - scale - amplification factor for the step that needs to be taken on actively bounded variables
330 
331   Output Parameter:
332 . S - step direction to be modified
333 
334   Level: developer
335 @*/
336 PetscErrorCode TaoBoundStep(Vec X, Vec XL, Vec XU, IS active_lower, IS active_upper, IS active_fixed, PetscReal scale, Vec S)
337 {
338 
339   Vec                          step_lower, step_upper, step_fixed;
340   Vec                          x_lower, x_upper;
341   Vec                          bound_lower, bound_upper;
342 
343   PetscFunctionBegin;
344   /* Adjust step for variables at the estimated lower bound */
345   if (active_lower) {
346     PetscCall(VecGetSubVector(S, active_lower, &step_lower));
347     PetscCall(VecGetSubVector(X, active_lower, &x_lower));
348     PetscCall(VecGetSubVector(XL, active_lower, &bound_lower));
349     PetscCall(VecCopy(bound_lower, step_lower));
350     PetscCall(VecAXPY(step_lower, -1.0, x_lower));
351     PetscCall(VecScale(step_lower, scale));
352     PetscCall(VecRestoreSubVector(S, active_lower, &step_lower));
353     PetscCall(VecRestoreSubVector(X, active_lower, &x_lower));
354     PetscCall(VecRestoreSubVector(XL, active_lower, &bound_lower));
355   }
356 
357   /* Adjust step for the variables at the estimated upper bound */
358   if (active_upper) {
359     PetscCall(VecGetSubVector(S, active_upper, &step_upper));
360     PetscCall(VecGetSubVector(X, active_upper, &x_upper));
361     PetscCall(VecGetSubVector(XU, active_upper, &bound_upper));
362     PetscCall(VecCopy(bound_upper, step_upper));
363     PetscCall(VecAXPY(step_upper, -1.0, x_upper));
364     PetscCall(VecScale(step_upper, scale));
365     PetscCall(VecRestoreSubVector(S, active_upper, &step_upper));
366     PetscCall(VecRestoreSubVector(X, active_upper, &x_upper));
367     PetscCall(VecRestoreSubVector(XU, active_upper, &bound_upper));
368   }
369 
370   /* Zero out step for fixed variables */
371   if (active_fixed) {
372     PetscCall(VecGetSubVector(S, active_fixed, &step_fixed));
373     PetscCall(VecSet(step_fixed, 0.0));
374     PetscCall(VecRestoreSubVector(S, active_fixed, &step_fixed));
375   }
376   PetscFunctionReturn(0);
377 }
378 
379 /*@C
380   TaoBoundSolution - Ensures that the solution vector is snapped into the bounds within a given tolerance.
381 
382   Collective on Vec
383 
384   Input Parameters:
385 + X - solution vector
386 . XL - lower bound vector
387 . XU - upper bound vector
388 - bound_tol - absolute tolerance in enforcing the bound
389 
390   Output Parameters:
391 + nDiff - total number of vector entries that have been bounded
392 - Xout - modified solution vector satisfying bounds to bound_tol
393 
394   Level: developer
395 
396 .seealso: `TAOBNCG`, `TAOBNTL`, `TAOBNTR`
397 @*/
398 PetscErrorCode TaoBoundSolution(Vec X, Vec XL, Vec XU, PetscReal bound_tol, PetscInt *nDiff, Vec Xout)
399 {
400   PetscInt          i,n,low,high,nDiff_loc=0;
401   PetscScalar       *xout;
402   const PetscScalar *x,*xl,*xu;
403 
404   PetscFunctionBegin;
405   PetscValidHeaderSpecific(X,VEC_CLASSID,1);
406   if (XL) PetscValidHeaderSpecific(XL,VEC_CLASSID,2);
407   if (XU) PetscValidHeaderSpecific(XU,VEC_CLASSID,3);
408   PetscValidHeaderSpecific(Xout,VEC_CLASSID,6);
409   if (!XL && !XU) {
410     PetscCall(VecCopy(X,Xout));
411     *nDiff = 0.0;
412     PetscFunctionReturn(0);
413   }
414   PetscCheckSameType(X,1,XL,2);
415   PetscCheckSameType(X,1,XU,3);
416   PetscCheckSameType(X,1,Xout,6);
417   PetscCheckSameComm(X,1,XL,2);
418   PetscCheckSameComm(X,1,XU,3);
419   PetscCheckSameComm(X,1,Xout,6);
420   VecCheckSameSize(X,1,XL,2);
421   VecCheckSameSize(X,1,XU,3);
422   VecCheckSameSize(X,1,Xout,4);
423 
424   PetscCall(VecGetOwnershipRange(X,&low,&high));
425   PetscCall(VecGetLocalSize(X,&n));
426   if (n>0) {
427     PetscCall(VecGetArrayRead(X, &x));
428     PetscCall(VecGetArrayRead(XL, &xl));
429     PetscCall(VecGetArrayRead(XU, &xu));
430     PetscCall(VecGetArray(Xout, &xout));
431 
432     for (i=0;i<n;++i) {
433       if (xl[i] > PETSC_NINFINITY && x[i] <= xl[i] + bound_tol) {
434         xout[i] = xl[i]; ++nDiff_loc;
435       } else if (xu[i] < PETSC_INFINITY && x[i] >= xu[i] - bound_tol) {
436         xout[i] = xu[i]; ++nDiff_loc;
437       }
438     }
439 
440     PetscCall(VecRestoreArrayRead(X, &x));
441     PetscCall(VecRestoreArrayRead(XL, &xl));
442     PetscCall(VecRestoreArrayRead(XU, &xu));
443     PetscCall(VecRestoreArray(Xout, &xout));
444   }
445   PetscCall(MPIU_Allreduce(&nDiff_loc, nDiff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)X)));
446   PetscFunctionReturn(0);
447 }
448