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