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