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