xref: /petsc/src/snes/impls/vi/vi.c (revision c3b5f7ba6bc5ce25a01a67bb37ba5d34b02bbbd7)
1 #include <petsc/private/snesimpl.h>  /*I "petscsnes.h" I*/
2 #include <petscdm.h>
3 
4 /*@C
5    SNESVISetComputeVariableBounds - Sets a function that is called to compute the variable bounds
6 
7    Input parameter:
8 +  snes - the SNES context
9 -  compute - computes the bounds
10 
11    Level: advanced
12 
13 .seealso:   SNESVISetVariableBounds()
14 
15 @*/
16 PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
17 {
18   PetscErrorCode (*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec));
19 
20   PetscFunctionBegin;
21   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
22   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",&f));
23   if (f) PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec)),(snes,compute));
24   else PetscCall(SNESVISetComputeVariableBounds_VI(snes,compute));
25   PetscFunctionReturn(0);
26 }
27 
28 PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
29 {
30   PetscFunctionBegin;
31   snes->ops->computevariablebounds = compute;
32   PetscFunctionReturn(0);
33 }
34 
35 /* --------------------------------------------------------------------------------------------------------*/
36 
37 PetscErrorCode  SNESVIMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
38 {
39   Vec            X, F, Finactive;
40   IS             isactive;
41   PetscViewer    viewer = (PetscViewer) dummy;
42 
43   PetscFunctionBegin;
44   PetscCall(SNESGetFunction(snes,&F,NULL,NULL));
45   PetscCall(SNESGetSolution(snes,&X));
46   PetscCall(SNESVIGetActiveSetIS(snes,X,F,&isactive));
47   PetscCall(VecDuplicate(F,&Finactive));
48   PetscCall(VecCopy(F,Finactive));
49   PetscCall(VecISSet(Finactive,isactive,0.0));
50   PetscCall(ISDestroy(&isactive));
51   PetscCall(VecView(Finactive,viewer));
52   PetscCall(VecDestroy(&Finactive));
53   PetscFunctionReturn(0);
54 }
55 
56 PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
57 {
58   PetscViewer       viewer = (PetscViewer) dummy;
59   const PetscScalar *x,*xl,*xu,*f;
60   PetscInt          i,n,act[2] = {0,0},fact[2],N;
61   /* Number of components that actually hit the bounds (c.f. active variables) */
62   PetscInt          act_bound[2] = {0,0},fact_bound[2];
63   PetscReal         rnorm,fnorm,zerotolerance = snes->vizerotolerance;
64   double            tmp;
65 
66   PetscFunctionBegin;
67   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
68   PetscCall(VecGetLocalSize(snes->vec_sol,&n));
69   PetscCall(VecGetSize(snes->vec_sol,&N));
70   PetscCall(VecGetArrayRead(snes->xl,&xl));
71   PetscCall(VecGetArrayRead(snes->xu,&xu));
72   PetscCall(VecGetArrayRead(snes->vec_sol,&x));
73   PetscCall(VecGetArrayRead(snes->vec_func,&f));
74 
75   rnorm = 0.0;
76   for (i=0; i<n; i++) {
77     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
78     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++;
79     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++;
80     else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Can never get here");
81   }
82 
83   for (i=0; i<n; i++) {
84     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++;
85     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++;
86   }
87   PetscCall(VecRestoreArrayRead(snes->vec_func,&f));
88   PetscCall(VecRestoreArrayRead(snes->xl,&xl));
89   PetscCall(VecRestoreArrayRead(snes->xu,&xu));
90   PetscCall(VecRestoreArrayRead(snes->vec_sol,&x));
91   PetscCall(MPIU_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes)));
92   PetscCall(MPIU_Allreduce(act,fact,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes)));
93   PetscCall(MPIU_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes)));
94   fnorm = PetscSqrtReal(fnorm);
95 
96   PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel));
97   if (snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds);
98   else tmp = 0.0;
99   PetscCall(PetscViewerASCIIPrintf(viewer,"%3D SNES VI Function norm %g Active lower constraints %D/%D upper constraints %D/%D Percent of total %g Percent of bounded %g\n",its,(double)fnorm,fact[0],fact_bound[0],fact[1],fact_bound[1],((double)(fact[0]+fact[1]))/((double)N),tmp));
100 
101   PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel));
102   PetscFunctionReturn(0);
103 }
104 
105 /*
106      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
107     || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that
108     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
109     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
110 */
111 PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
112 {
113   PetscReal      a1;
114   PetscBool      hastranspose;
115 
116   PetscFunctionBegin;
117   *ismin = PETSC_FALSE;
118   PetscCall(MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose));
119   if (hastranspose) {
120     /* Compute || J^T F|| */
121     PetscCall(MatMultTranspose(A,F,W));
122     PetscCall(VecNorm(W,NORM_2,&a1));
123     PetscCall(PetscInfo(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm)));
124     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
125   } else {
126     Vec         work;
127     PetscScalar result;
128     PetscReal   wnorm;
129 
130     PetscCall(VecSetRandom(W,NULL));
131     PetscCall(VecNorm(W,NORM_2,&wnorm));
132     PetscCall(VecDuplicate(W,&work));
133     PetscCall(MatMult(A,W,work));
134     PetscCall(VecDot(F,work,&result));
135     PetscCall(VecDestroy(&work));
136     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
137     PetscCall(PetscInfo(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1));
138     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
139   }
140   PetscFunctionReturn(0);
141 }
142 
143 /*
144      Checks if J^T(F - J*X) = 0
145 */
146 PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
147 {
148   PetscReal      a1,a2;
149   PetscBool      hastranspose;
150 
151   PetscFunctionBegin;
152   PetscCall(MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose));
153   if (hastranspose) {
154     PetscCall(MatMult(A,X,W1));
155     PetscCall(VecAXPY(W1,-1.0,F));
156 
157     /* Compute || J^T W|| */
158     PetscCall(MatMultTranspose(A,W1,W2));
159     PetscCall(VecNorm(W1,NORM_2,&a1));
160     PetscCall(VecNorm(W2,NORM_2,&a2));
161     if (a1 != 0.0) {
162       PetscCall(PetscInfo(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1)));
163     }
164   }
165   PetscFunctionReturn(0);
166 }
167 
168 /*
169   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
170 
171   Notes:
172   The convergence criterion currently implemented is
173   merit < abstol
174   merit < rtol*merit_initial
175 */
176 PetscErrorCode SNESConvergedDefault_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
177 {
178   PetscFunctionBegin;
179   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
180   PetscValidPointer(reason,6);
181 
182   *reason = SNES_CONVERGED_ITERATING;
183 
184   if (!it) {
185     /* set parameter for default relative tolerance convergence test */
186     snes->ttol = fnorm*snes->rtol;
187   }
188   if (fnorm != fnorm) {
189     PetscCall(PetscInfo(snes,"Failed to converged, function norm is NaN\n"));
190     *reason = SNES_DIVERGED_FNORM_NAN;
191   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
192     PetscCall(PetscInfo(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol));
193     *reason = SNES_CONVERGED_FNORM_ABS;
194   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
195     PetscCall(PetscInfo(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs));
196     *reason = SNES_DIVERGED_FUNCTION_COUNT;
197   }
198 
199   if (it && !*reason) {
200     if (fnorm < snes->ttol) {
201       PetscCall(PetscInfo(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol));
202       *reason = SNES_CONVERGED_FNORM_RELATIVE;
203     }
204   }
205   PetscFunctionReturn(0);
206 }
207 
208 /* -------------------------------------------------------------------------- */
209 /*
210    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
211 
212    Input Parameters:
213 .  SNES - nonlinear solver context
214 
215    Output Parameters:
216 .  X - Bound projected X
217 
218 */
219 
220 PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
221 {
222   const PetscScalar *xl,*xu;
223   PetscScalar       *x;
224   PetscInt          i,n;
225 
226   PetscFunctionBegin;
227   PetscCall(VecGetLocalSize(X,&n));
228   PetscCall(VecGetArray(X,&x));
229   PetscCall(VecGetArrayRead(snes->xl,&xl));
230   PetscCall(VecGetArrayRead(snes->xu,&xu));
231 
232   for (i = 0; i<n; i++) {
233     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
234     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
235   }
236   PetscCall(VecRestoreArray(X,&x));
237   PetscCall(VecRestoreArrayRead(snes->xl,&xl));
238   PetscCall(VecRestoreArrayRead(snes->xu,&xu));
239   PetscFunctionReturn(0);
240 }
241 
242 /*
243    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
244 
245    Input parameter:
246 .  snes - the SNES context
247 .  X    - the snes solution vector
248 .  F    - the nonlinear function vector
249 
250    Output parameter:
251 .  ISact - active set index set
252  */
253 PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS *ISact)
254 {
255   Vec               Xl=snes->xl,Xu=snes->xu;
256   const PetscScalar *x,*f,*xl,*xu;
257   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
258   PetscReal         zerotolerance = snes->vizerotolerance;
259 
260   PetscFunctionBegin;
261   PetscCall(VecGetLocalSize(X,&nlocal));
262   PetscCall(VecGetOwnershipRange(X,&ilow,&ihigh));
263   PetscCall(VecGetArrayRead(X,&x));
264   PetscCall(VecGetArrayRead(Xl,&xl));
265   PetscCall(VecGetArrayRead(Xu,&xu));
266   PetscCall(VecGetArrayRead(F,&f));
267   /* Compute active set size */
268   for (i=0; i < nlocal;i++) {
269     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) nloc_isact++;
270   }
271 
272   PetscCall(PetscMalloc1(nloc_isact,&idx_act));
273 
274   /* Set active set indices */
275   for (i=0; i < nlocal; i++) {
276     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) idx_act[i1++] = ilow+i;
277   }
278 
279   /* Create active set IS */
280   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes),nloc_isact,idx_act,PETSC_OWN_POINTER,ISact));
281 
282   PetscCall(VecRestoreArrayRead(X,&x));
283   PetscCall(VecRestoreArrayRead(Xl,&xl));
284   PetscCall(VecRestoreArrayRead(Xu,&xu));
285   PetscCall(VecRestoreArrayRead(F,&f));
286   PetscFunctionReturn(0);
287 }
288 
289 PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact)
290 {
291   PetscInt       rstart,rend;
292 
293   PetscFunctionBegin;
294   PetscCall(SNESVIGetActiveSetIS(snes,X,F,ISact));
295   PetscCall(VecGetOwnershipRange(X,&rstart,&rend));
296   PetscCall(ISComplement(*ISact,rstart,rend,ISinact));
297   PetscFunctionReturn(0);
298 }
299 
300 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm)
301 {
302   const PetscScalar *x,*xl,*xu,*f;
303   PetscInt          i,n;
304   PetscReal         rnorm,zerotolerance = snes->vizerotolerance;
305 
306   PetscFunctionBegin;
307   PetscCall(VecGetLocalSize(X,&n));
308   PetscCall(VecGetArrayRead(snes->xl,&xl));
309   PetscCall(VecGetArrayRead(snes->xu,&xu));
310   PetscCall(VecGetArrayRead(X,&x));
311   PetscCall(VecGetArrayRead(F,&f));
312   rnorm = 0.0;
313   for (i=0; i<n; i++) {
314     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
315   }
316   PetscCall(VecRestoreArrayRead(F,&f));
317   PetscCall(VecRestoreArrayRead(snes->xl,&xl));
318   PetscCall(VecRestoreArrayRead(snes->xu,&xu));
319   PetscCall(VecRestoreArrayRead(X,&x));
320   PetscCall(MPIU_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes)));
321   *fnorm = PetscSqrtReal(*fnorm);
322   PetscFunctionReturn(0);
323 }
324 
325 PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
326 {
327   PetscFunctionBegin;
328   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
329   PetscFunctionReturn(0);
330 }
331 
332 /* -------------------------------------------------------------------------- */
333 /*
334    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
335    of the SNESVI nonlinear solver.
336 
337    Input Parameter:
338 .  snes - the SNES context
339 
340    Application Interface Routine: SNESSetUp()
341 
342    Notes:
343    For basic use of the SNES solvers, the user need not explicitly call
344    SNESSetUp(), since these actions will automatically occur during
345    the call to SNESSolve().
346  */
347 PetscErrorCode SNESSetUp_VI(SNES snes)
348 {
349   PetscInt       i_start[3],i_end[3];
350 
351   PetscFunctionBegin;
352   PetscCall(SNESSetWorkVecs(snes,1));
353   PetscCall(SNESSetUpMatrices(snes));
354 
355   if (!snes->ops->computevariablebounds && snes->dm) {
356     PetscBool flag;
357     PetscCall(DMHasVariableBounds(snes->dm, &flag));
358     if (flag) {
359       snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
360     }
361   }
362   if (!snes->usersetbounds) {
363     if (snes->ops->computevariablebounds) {
364       if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol,&snes->xl));
365       if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol,&snes->xu));
366       PetscCall((*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu));
367     } else if (!snes->xl && !snes->xu) {
368       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
369       PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
370       PetscCall(VecSet(snes->xl,PETSC_NINFINITY));
371       PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
372       PetscCall(VecSet(snes->xu,PETSC_INFINITY));
373     } else {
374       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
375       PetscCall(VecGetOwnershipRange(snes->vec_sol,i_start,i_end));
376       PetscCall(VecGetOwnershipRange(snes->xl,i_start+1,i_end+1));
377       PetscCall(VecGetOwnershipRange(snes->xu,i_start+2,i_end+2));
378       if ((i_start[0] != i_start[1]) || (i_start[0] != i_start[2]) || (i_end[0] != i_end[1]) || (i_end[0] != i_end[2]))
379         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Distribution of lower bound, upper bound and the solution vector should be identical across all the processors.");
380     }
381   }
382   PetscFunctionReturn(0);
383 }
384 /* -------------------------------------------------------------------------- */
385 PetscErrorCode SNESReset_VI(SNES snes)
386 {
387   PetscFunctionBegin;
388   PetscCall(VecDestroy(&snes->xl));
389   PetscCall(VecDestroy(&snes->xu));
390   snes->usersetbounds = PETSC_FALSE;
391   PetscFunctionReturn(0);
392 }
393 
394 /*
395    SNESDestroy_VI - Destroys the private SNES_VI context that was created
396    with SNESCreate_VI().
397 
398    Input Parameter:
399 .  snes - the SNES context
400 
401    Application Interface Routine: SNESDestroy()
402  */
403 PetscErrorCode SNESDestroy_VI(SNES snes)
404 {
405   PetscFunctionBegin;
406   PetscCall(PetscFree(snes->data));
407 
408   /* clear composed functions */
409   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSet_C",NULL));
410   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSetDefaultMonitor_C",NULL));
411   PetscFunctionReturn(0);
412 }
413 
414 /*@
415    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
416 
417    Input Parameters:
418 +  snes - the SNES context.
419 .  xl   - lower bound.
420 -  xu   - upper bound.
421 
422    Notes:
423    If this routine is not called then the lower and upper bounds are set to
424    PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp().
425 
426    Level: advanced
427 
428 @*/
429 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
430 {
431   PetscErrorCode (*f)(SNES,Vec,Vec);
432 
433   PetscFunctionBegin;
434   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
435   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
436   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
437   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",&f));
438   if (f) PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));
439   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
440   snes->usersetbounds = PETSC_TRUE;
441   PetscFunctionReturn(0);
442 }
443 
444 PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
445 {
446   const PetscScalar *xxl,*xxu;
447   PetscInt          i,n, cnt = 0;
448 
449   PetscFunctionBegin;
450   PetscCall(SNESGetFunction(snes,&snes->vec_func,NULL,NULL));
451   PetscCheck(snes->vec_func,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
452   {
453     PetscInt xlN,xuN,N;
454     PetscCall(VecGetSize(xl,&xlN));
455     PetscCall(VecGetSize(xu,&xuN));
456     PetscCall(VecGetSize(snes->vec_func,&N));
457     PetscCheck(xlN == N,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xlN,N);
458     PetscCheck(xuN == N,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xuN,N);
459   }
460   PetscCall(PetscObjectReference((PetscObject)xl));
461   PetscCall(PetscObjectReference((PetscObject)xu));
462   PetscCall(VecDestroy(&snes->xl));
463   PetscCall(VecDestroy(&snes->xu));
464   snes->xl = xl;
465   snes->xu = xu;
466   PetscCall(VecGetLocalSize(xl,&n));
467   PetscCall(VecGetArrayRead(xl,&xxl));
468   PetscCall(VecGetArrayRead(xu,&xxu));
469   for (i=0; i<n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
470 
471   PetscCall(MPIU_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes)));
472   PetscCall(VecRestoreArrayRead(xl,&xxl));
473   PetscCall(VecRestoreArrayRead(xu,&xxu));
474   PetscFunctionReturn(0);
475 }
476 
477 PetscErrorCode SNESSetFromOptions_VI(PetscOptionItems *PetscOptionsObject,SNES snes)
478 {
479   PetscBool      flg = PETSC_FALSE;
480 
481   PetscFunctionBegin;
482   PetscCall(PetscOptionsHead(PetscOptionsObject,"SNES VI options"));
483   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance","Tolerance for considering x[] value to be on a bound","None",snes->vizerotolerance,&snes->vizerotolerance,NULL));
484   PetscCall(PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","SNESMonitorResidual",flg,&flg,NULL));
485   if (flg) {
486     PetscCall(SNESMonitorSet(snes,SNESMonitorVI,PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),NULL));
487   }
488   flg = PETSC_FALSE;
489   PetscCall(PetscOptionsBool("-snes_vi_monitor_residual","Monitor residual all non-active variables; using zero for active constraints","SNESMonitorVIResidual",flg,&flg,NULL));
490   if (flg) {
491     PetscCall(SNESMonitorSet(snes,SNESVIMonitorResidual,PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)),NULL));
492   }
493   PetscCall(PetscOptionsTail());
494   PetscFunctionReturn(0);
495 }
496