xref: /petsc/src/snes/impls/vi/vi.c (revision 1db2f0e57964f94f63a3f5cf6bb1e5befdb56f30)
1 
2 #include <../include/private/snesimpl.h> /*I "petscsnes.h" I*/
3 #include <../include/private/kspimpl.h>
4 #include <../include/private/matimpl.h>
5 #include <../include/private/dmimpl.h>
6 
7 #undef __FUNCT__
8 #define __FUNCT__ "SNESVISetComputeVariableBounds"
9 /*@C
10    SNESVISetComputeVariableBounds - Sets a function  that is called to compute the variable bounds
11 
12    Input parameter
13 +  snes - the SNES context
14 -  compute - computes the bounds
15 
16    Level: advanced
17 
18 @*/
19 PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
20 {
21   PetscErrorCode   ierr,(*f)(SNES,PetscErrorCode(*)(SNES,Vec,Vec));
22 
23   PetscFunctionBegin;
24   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
25   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
26   if (!f) {ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);}
27   ierr = PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode(*)(SNES,Vec,Vec)),(snes,compute));CHKERRQ(ierr);
28   PetscFunctionReturn(0);
29 }
30 
31 EXTERN_C_BEGIN
32 #undef __FUNCT__
33 #define __FUNCT__ "SNESVISetComputeVariableBounds_VI"
34 PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
35 {
36   PetscFunctionBegin;
37   snes->ops->computevariablebounds = compute;
38   PetscFunctionReturn(0);
39 }
40 EXTERN_C_END
41 
42 #undef __FUNCT__
43 #define __FUNCT__ "SNESVIComputeInactiveSetIS"
44 /*
45    SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables
46 
47    Input parameter
48 .  snes - the SNES context
49 .  X    - the snes solution vector
50 
51    Output parameter
52 .  ISact - active set index set
53 
54  */
55 PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS* inact)
56 {
57   PetscErrorCode    ierr;
58   const PetscScalar *x,*xl,*xu,*f;
59   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
60 
61   PetscFunctionBegin;
62   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
63   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
64   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
65   ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr);
66   ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr);
67   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
68   /* Compute inactive set size */
69   for (i=0; i < nlocal;i++) {
70     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++;
71   }
72 
73   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
74 
75   /* Set inactive set indices */
76   for(i=0; i < nlocal; i++) {
77     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;
78   }
79 
80    /* Create inactive set IS */
81   ierr = ISCreateGeneral(((PetscObject)upper)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr);
82 
83   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
84   ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr);
85   ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr);
86   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 /* --------------------------------------------------------------------------------------------------------*/
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "SNESMonitorVI"
94 PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
95 {
96   PetscErrorCode    ierr;
97   PetscViewer        viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
98   const PetscScalar  *x,*xl,*xu,*f;
99   PetscInt           i,n,act[2] = {0,0},fact[2],N;
100   /* Number of components that actually hit the bounds (c.f. active variables) */
101   PetscInt           act_bound[2] = {0,0},fact_bound[2];
102   PetscReal          rnorm,fnorm;
103 
104   PetscFunctionBegin;
105   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
106   ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr);
107   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
108   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
109   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
110   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
111 
112   rnorm = 0.0;
113   for (i=0; i<n; i++) {
114     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]);
115     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++;
116     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++;
117     else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Can never get here");
118   }
119 
120   for (i=0; i<n; i++) {
121     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++;
122     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++;
123   }
124   ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
125   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
126   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
127   ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
128   ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
129   ierr = MPI_Allreduce(act,fact,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
130   ierr = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
131   fnorm = PetscSqrtReal(fnorm);
132 
133   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
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),((double)(fact[0]+fact[1]))/((double)snes->ntruebounds));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,PETSC_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   SNESDefaultConverged_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__ "SNESDefaultConverged_VI"
219 PetscErrorCode SNESDefaultConverged_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 = PetscMalloc(nloc_isact*sizeof(PetscInt),&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(((PetscObject)snes)->comm,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 
347   PetscFunctionBegin;
348   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
349   ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr);
350   PetscFunctionReturn(0);
351 }
352 
353 #undef __FUNCT__
354 #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
355 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X,PetscReal *fnorm)
356 {
357   PetscErrorCode    ierr;
358   const PetscScalar *x,*xl,*xu,*f;
359   PetscInt          i,n;
360   PetscReal         rnorm;
361 
362   PetscFunctionBegin;
363   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
364   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
365   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
366   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
367   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
368   rnorm = 0.0;
369   for (i=0; i<n; i++) {
370     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]);
371   }
372   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
373   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
374   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
375   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
376   ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
377   *fnorm = PetscSqrtReal(*fnorm);
378   PetscFunctionReturn(0);
379 }
380 
381 #undef __FUNCT__
382 #define __FUNCT__ "SNESVIDMComputeVariableBounds"
383 PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
384 {
385   PetscErrorCode     ierr;
386   PetscFunctionBegin;
387   ierr = DMComputeVariableBounds(snes->dm, xl, xu); CHKERRQ(ierr);
388   PetscFunctionReturn(0);
389 }
390 
391 
392 /* -------------------------------------------------------------------------- */
393 /*
394    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
395    of the SNESVI nonlinear solver.
396 
397    Input Parameter:
398 .  snes - the SNES context
399 .  x - the solution vector
400 
401    Application Interface Routine: SNESSetUp()
402 
403    Notes:
404    For basic use of the SNES solvers, the user need not explicitly call
405    SNESSetUp(), since these actions will automatically occur during
406    the call to SNESSolve().
407  */
408 #undef __FUNCT__
409 #define __FUNCT__ "SNESSetUp_VI"
410 PetscErrorCode SNESSetUp_VI(SNES snes)
411 {
412   PetscErrorCode ierr;
413   PetscInt       i_start[3],i_end[3];
414 
415   PetscFunctionBegin;
416 
417   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
418   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
419 
420   if(!snes->ops->computevariablebounds && snes->dm) {
421     snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
422   }
423   if (snes->ops->computevariablebounds) {
424     if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
425     if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
426     ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
427   } else if (!snes->xl && !snes->xu) {
428     /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
429     ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr);
430     ierr = VecSet(snes->xl,SNES_VI_NINF);CHKERRQ(ierr);
431     ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr);
432     ierr = VecSet(snes->xu,SNES_VI_INF);CHKERRQ(ierr);
433   } else {
434     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
435     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
436     ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr);
437     ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr);
438     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]))
439       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.");
440   }
441 
442   PetscFunctionReturn(0);
443 }
444 /* -------------------------------------------------------------------------- */
445 #undef __FUNCT__
446 #define __FUNCT__ "SNESReset_VI"
447 PetscErrorCode SNESReset_VI(SNES snes)
448 {
449   PetscErrorCode ierr;
450 
451   PetscFunctionBegin;
452   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
453   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
454   PetscFunctionReturn(0);
455 }
456 
457 /*
458    SNESDestroy_VI - Destroys the private SNES_VI context that was created
459    with SNESCreate_VI().
460 
461    Input Parameter:
462 .  snes - the SNES context
463 
464    Application Interface Routine: SNESDestroy()
465  */
466 #undef __FUNCT__
467 #define __FUNCT__ "SNESDestroy_VI"
468 PetscErrorCode SNESDestroy_VI(SNES snes)
469 {
470   PetscErrorCode ierr;
471 
472   PetscFunctionBegin;
473   ierr = PetscFree(snes->data);CHKERRQ(ierr);
474 
475   /* clear composed functions */
476   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
477   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
478   PetscFunctionReturn(0);
479 }
480 
481 /* -------------------------------------------------------------------------- */
482 
483 /*
484 
485       These line searches are common for all the VI solvers
486 */
487 extern PetscErrorCode SNESSolve_VISS(SNES);
488 
489 #undef __FUNCT__
490 #define __FUNCT__ "SNESLineSearchNo_VI"
491 
492 /*
493   This routine does not actually do a line search but it takes a full newton
494   step while ensuring that the new iterates remain within the constraints.
495 
496 */
497 PetscErrorCode SNESLineSearchNo_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
498 {
499   PetscErrorCode ierr;
500   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
501 
502   PetscFunctionBegin;
503   *flag = PETSC_TRUE;
504   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
505   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
506   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
507   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
508   if (snes->ops->postcheckstep) {
509    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
510   }
511   if (changed_y) {
512     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
513     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
514   }
515   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
516   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
517   if (!snes->domainerror) {
518     if (snes->ops->solve != SNESSolve_VISS) {
519        ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
520     } else {
521       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
522     }
523     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
524   }
525   if (snes->ls_monitor) {
526     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
527     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
528     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
529   }
530   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
531   PetscFunctionReturn(0);
532 }
533 
534 /* -------------------------------------------------------------------------- */
535 #undef __FUNCT__
536 #define __FUNCT__ "SNESLineSearchNoNorms_VI"
537 
538 /*
539   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
540 */
541 PetscErrorCode SNESLineSearchNoNorms_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
542 {
543   PetscErrorCode ierr;
544   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
545 
546   PetscFunctionBegin;
547   *flag = PETSC_TRUE;
548   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
549   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
550   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
551   if (snes->ops->postcheckstep) {
552    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
553   }
554   if (changed_y) {
555     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
556     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
557   }
558 
559   /* don't evaluate function the last time through */
560   if (snes->iter < snes->max_its-1) {
561     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
562   }
563   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
564   PetscFunctionReturn(0);
565 }
566 
567 /* -------------------------------------------------------------------------- */
568 
569 #undef __FUNCT__
570 #define __FUNCT__ "SNESLineSearchCubic_VI"
571 /*
572   This routine implements a cubic line search while doing a projection on the variable bounds
573 */
574 PetscErrorCode SNESLineSearchCubic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
575 {
576   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
577   PetscReal      minlambda,lambda,lambdatemp;
578 #if defined(PETSC_USE_COMPLEX)
579   PetscScalar    cinitslope;
580 #endif
581   PetscErrorCode ierr;
582   PetscInt       count;
583   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
584   MPI_Comm       comm;
585 
586   PetscFunctionBegin;
587   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
588   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
589   *flag   = PETSC_TRUE;
590 
591   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
592   if (*ynorm == 0.0) {
593     if (snes->ls_monitor) {
594       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
595       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
596       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
597     }
598     *gnorm = fnorm;
599     ierr   = VecCopy(x,w);CHKERRQ(ierr);
600     ierr   = VecCopy(f,g);CHKERRQ(ierr);
601     *flag  = PETSC_FALSE;
602     goto theend1;
603   }
604   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
605     if (snes->ls_monitor) {
606       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
607       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Scaling step by %g old ynorm %g\n",(double)(snes->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr);
608       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
609     }
610     ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
611     *ynorm = snes->maxstep;
612   }
613   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
614   minlambda = snes->steptol/rellength;
615   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
616 #if defined(PETSC_USE_COMPLEX)
617   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
618   initslope = PetscRealPart(cinitslope);
619 #else
620   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
621 #endif
622   if (initslope > 0.0)  initslope = -initslope;
623   if (initslope == 0.0) initslope = -1.0;
624 
625   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
626   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
627   if (snes->nfuncs >= snes->max_funcs) {
628     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
629     *flag = PETSC_FALSE;
630     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
631     goto theend1;
632   }
633   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
634   if (snes->ops->solve != SNESSolve_VISS) {
635     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
636   } else {
637     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
638   }
639   if (snes->domainerror) {
640     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
641     PetscFunctionReturn(0);
642   }
643   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
644   ierr = PetscInfo4(snes,"Initial fnorm %g gnorm %g alpha %g initslope %g\n",(double)fnorm,(double)*gnorm,(double)snes->ls_alpha,(double)initslope);CHKERRQ(ierr);
645   if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* Sufficient reduction */
646     if (snes->ls_monitor) {
647       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
648       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr);
649       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
650     }
651     goto theend1;
652   }
653 
654   /* Fit points with quadratic */
655   lambda     = 1.0;
656   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
657   lambdaprev = lambda;
658   gnormprev  = *gnorm;
659   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
660   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
661   else                         lambda = lambdatemp;
662 
663   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
664   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
665   if (snes->nfuncs >= snes->max_funcs) {
666     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
667     *flag = PETSC_FALSE;
668     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
669     goto theend1;
670   }
671   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
672   if (snes->ops->solve != SNESSolve_VISS) {
673     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
674   } else {
675     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
676   }
677   if (snes->domainerror) {
678     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
679     PetscFunctionReturn(0);
680   }
681   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
682   if (snes->ls_monitor) {
683     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
684     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: gnorm after quadratic fit %g\n",(double)(*gnorm));CHKERRQ(ierr);
685     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
686   }
687   if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* sufficient reduction */
688     if (snes->ls_monitor) {
689       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
690       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
691       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
692     }
693     goto theend1;
694   }
695 
696   /* Fit points with cubic */
697   count = 1;
698   while (PETSC_TRUE) {
699     if (lambda <= minlambda) {
700       if (snes->ls_monitor) {
701         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
702  	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
703 	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)minlambda,(double)lambda,(double)initslope);CHKERRQ(ierr);
704         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
705       }
706       *flag = PETSC_FALSE;
707       break;
708     }
709     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
710     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
711     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
712     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
713     d  = b*b - 3*a*initslope;
714     if (d < 0.0) d = 0.0;
715     if (a == 0.0) {
716       lambdatemp = -initslope/(2.0*b);
717     } else {
718       lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
719     }
720     lambdaprev = lambda;
721     gnormprev  = *gnorm;
722     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
723     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
724     else                         lambda     = lambdatemp;
725     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
726     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
727     if (snes->nfuncs >= snes->max_funcs) {
728       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
729       ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)lambda,(double)initslope);CHKERRQ(ierr);
730       *flag = PETSC_FALSE;
731       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
732       break;
733     }
734     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
735     if (snes->ops->solve != SNESSolve_VISS) {
736       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
737     } else {
738       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
739     }
740     if (snes->domainerror) {
741       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
742       PetscFunctionReturn(0);
743     }
744     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
745     if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* is reduction enough? */
746       if (snes->ls_monitor) {
747 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %g lambda=%18.16e\n",(double)(*gnorm),(double)lambda);CHKERRQ(ierr);
748       }
749       break;
750     } else {
751       if (snes->ls_monitor) {
752         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %g lambda=%18.16e\n",(double)(*gnorm),(double)lambda);CHKERRQ(ierr);
753       }
754     }
755     count++;
756   }
757   theend1:
758   /* Optional user-defined check for line search step validity */
759   if (snes->ops->postcheckstep && *flag) {
760     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
761     if (changed_y) {
762       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
763       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
764     }
765     if (changed_y || changed_w) { /* recompute the function if the step has changed */
766       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
767       if (snes->ops->solve != SNESSolve_VISS) {
768         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
769       } else {
770         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
771       }
772       if (snes->domainerror) {
773         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
774         PetscFunctionReturn(0);
775       }
776       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
777       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
778       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
779     }
780   }
781   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
782   PetscFunctionReturn(0);
783 }
784 
785 #undef __FUNCT__
786 #define __FUNCT__ "SNESLineSearchQuadratic_VI"
787 /*
788   This routine does a quadratic line search while keeping the iterates within the variable bounds
789 */
790 PetscErrorCode SNESLineSearchQuadratic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
791 {
792   /*
793      Note that for line search purposes we work with with the related
794      minimization problem:
795         min  z(x):  R^n -> R,
796      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
797    */
798   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
799 #if defined(PETSC_USE_COMPLEX)
800   PetscScalar    cinitslope;
801 #endif
802   PetscErrorCode ierr;
803   PetscInt       count;
804   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
805 
806   PetscFunctionBegin;
807   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
808   *flag   = PETSC_TRUE;
809 
810   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
811   if (*ynorm == 0.0) {
812     if (snes->ls_monitor) {
813       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
814       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
815       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
816     }
817     *gnorm = fnorm;
818     ierr   = VecCopy(x,w);CHKERRQ(ierr);
819     ierr   = VecCopy(f,g);CHKERRQ(ierr);
820     *flag  = PETSC_FALSE;
821     goto theend2;
822   }
823   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
824     ierr   = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
825     *ynorm = snes->maxstep;
826   }
827   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
828   minlambda = snes->steptol/rellength;
829   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
830 #if defined(PETSC_USE_COMPLEX)
831   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
832   initslope = PetscRealPart(cinitslope);
833 #else
834   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
835 #endif
836   if (initslope > 0.0)  initslope = -initslope;
837   if (initslope == 0.0) initslope = -1.0;
838   ierr = PetscInfo1(snes,"Initslope %g \n",(double)initslope);CHKERRQ(ierr);
839 
840   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
841   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
842   if (snes->nfuncs >= snes->max_funcs) {
843     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
844     *flag = PETSC_FALSE;
845     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
846     goto theend2;
847   }
848   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
849   if (snes->ops->solve != SNESSolve_VISS) {
850     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
851   } else {
852     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
853   }
854   if (snes->domainerror) {
855     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
856     PetscFunctionReturn(0);
857   }
858   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
859   if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* Sufficient reduction */
860     if (snes->ls_monitor) {
861       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
862       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr);
863       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
864     }
865     goto theend2;
866   }
867 
868   /* Fit points with quadratic */
869   lambda = 1.0;
870   count = 1;
871   while (PETSC_TRUE) {
872     if (lambda <= minlambda) { /* bad luck; use full step */
873       if (snes->ls_monitor) {
874         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
875         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
876         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)lambda,(double)initslope);CHKERRQ(ierr);
877         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
878       }
879       ierr = VecCopy(x,w);CHKERRQ(ierr);
880       *flag = PETSC_FALSE;
881       break;
882     }
883     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
884     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
885     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
886     else                         lambda     = lambdatemp;
887 
888     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
889     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
890     if (snes->nfuncs >= snes->max_funcs) {
891       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
892       ierr  = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)lambda,(double)initslope);CHKERRQ(ierr);
893       *flag = PETSC_FALSE;
894       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
895       break;
896     }
897     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
898     if (snes->domainerror) {
899       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
900       PetscFunctionReturn(0);
901     }
902     if (snes->ops->solve != SNESSolve_VISS) {
903       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
904     } else {
905       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
906     }
907     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
908     if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* sufficient reduction */
909       if (snes->ls_monitor) {
910         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
911         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Quadratically determined step, lambda=%g\n",(double)lambda);CHKERRQ(ierr);
912         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
913       }
914       break;
915     }
916     count++;
917   }
918   theend2:
919   /* Optional user-defined check for line search step validity */
920   if (snes->ops->postcheckstep) {
921     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
922     if (changed_y) {
923       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
924       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
925     }
926     if (changed_y || changed_w) { /* recompute the function if the step has changed */
927       ierr = SNESComputeFunction(snes,w,g);
928       if (snes->domainerror) {
929         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
930         PetscFunctionReturn(0);
931       }
932       if (snes->ops->solve != SNESSolve_VISS) {
933         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
934       } else {
935         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
936       }
937 
938       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
939       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
940       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
941     }
942   }
943   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
944   PetscFunctionReturn(0);
945 }
946 
947 
948 #undef __FUNCT__
949 #define __FUNCT__ "SNESVISetVariableBounds"
950 /*@
951    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
952 
953    Input Parameters:
954 .  snes - the SNES context.
955 .  xl   - lower bound.
956 .  xu   - upper bound.
957 
958    Notes:
959    If this routine is not called then the lower and upper bounds are set to
960    SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp().
961 
962    Level: advanced
963 
964 @*/
965 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
966 {
967   PetscErrorCode   ierr,(*f)(SNES,Vec,Vec);
968 
969   PetscFunctionBegin;
970   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
971   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
972   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
973   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
974   if (!f) {ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);}
975   ierr = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr);
976   PetscFunctionReturn(0);
977 }
978 
979 EXTERN_C_BEGIN
980 #undef __FUNCT__
981 #define __FUNCT__ "SNESVISetVariableBounds_VI"
982 PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
983 {
984   PetscErrorCode    ierr;
985   const PetscScalar *xxl,*xxu;
986   PetscInt          i,n, cnt = 0;
987 
988   PetscFunctionBegin;
989   ierr = SNESGetFunction(snes,&snes->vec_func,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
990   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
991   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);
992   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);
993   ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);
994   ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
995   ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
996   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
997   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
998   snes->xl = xl;
999   snes->xu = xu;
1000   ierr = VecGetLocalSize(xl,&n);CHKERRQ(ierr);
1001   ierr = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr);
1002   ierr = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr);
1003   for (i=0; i<n; i++) {
1004     cnt += ((xxl[i] != SNES_VI_NINF) || (xxu[i] != SNES_VI_INF));
1005   }
1006   ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
1007   ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr);
1008   ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr);
1009   PetscFunctionReturn(0);
1010 }
1011 EXTERN_C_END
1012 
1013 EXTERN_C_BEGIN
1014 #undef __FUNCT__
1015 #define __FUNCT__ "SNESLineSearchSetType_VI"
1016 PetscErrorCode  SNESLineSearchSetType_VI(SNES snes, SNESLineSearchType type)
1017 {
1018   PetscErrorCode ierr;
1019   PetscFunctionBegin;
1020 
1021   switch (type) {
1022   case SNES_LS_BASIC:
1023     ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
1024     break;
1025   case SNES_LS_BASIC_NONORMS:
1026     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
1027     break;
1028   case SNES_LS_QUADRATIC:
1029     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
1030     break;
1031   case SNES_LS_CUBIC:
1032     ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
1033     break;
1034   default:
1035     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
1036     break;
1037   }
1038   snes->ls_type = type;
1039   PetscFunctionReturn(0);
1040 }
1041 EXTERN_C_END
1042 
1043 #undef __FUNCT__
1044 #define __FUNCT__ "SNESSetFromOptions_VI"
1045 PetscErrorCode SNESSetFromOptions_VI(SNES snes)
1046 {
1047   PetscErrorCode  ierr;
1048   PetscBool       flg;
1049 
1050   PetscFunctionBegin;
1051   ierr = PetscOptionsHead("SNES VI options");CHKERRQ(ierr);
1052   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
1053   if (flg) {
1054     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
1055   }
1056   ierr = PetscOptionsTail();CHKERRQ(ierr);
1057   PetscFunctionReturn(0);
1058 }
1059