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