xref: /petsc/src/snes/impls/vi/vi.c (revision 3e336fdf4e375df90ef4b65b0f00c0c010d980fd)
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;
22 
23   PetscFunctionBegin;
24   ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);
25   snes->ops->computevariablebounds = compute;
26   PetscFunctionReturn(0);
27 }
28 
29 
30 #undef __FUNCT__
31 #define __FUNCT__ "SNESVIComputeInactiveSetIS"
32 /*
33    SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables
34 
35    Input parameter
36 .  snes - the SNES context
37 .  X    - the snes solution vector
38 
39    Output parameter
40 .  ISact - active set index set
41 
42  */
43 PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS* inact)
44 {
45   PetscErrorCode    ierr;
46   const PetscScalar *x,*xl,*xu,*f;
47   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
48 
49   PetscFunctionBegin;
50   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
51   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
52   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
53   ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr);
54   ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr);
55   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
56   /* Compute inactive set size */
57   for (i=0; i < nlocal;i++) {
58     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++;
59   }
60 
61   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
62 
63   /* Set inactive set indices */
64   for(i=0; i < nlocal; i++) {
65     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;
66   }
67 
68    /* Create inactive set IS */
69   ierr = ISCreateGeneral(((PetscObject)upper)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr);
70 
71   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
72   ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr);
73   ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr);
74   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
75   PetscFunctionReturn(0);
76 }
77 
78 /* --------------------------------------------------------------------------------------------------------*/
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "SNESMonitorVI"
82 PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
83 {
84   PetscErrorCode    ierr;
85   PetscViewer        viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
86   const PetscScalar  *x,*xl,*xu,*f;
87   PetscInt           i,n,act[2] = {0,0},fact[2],N;
88   /* Number of components that actually hit the bounds (c.f. active variables) */
89   PetscInt           act_bound[2] = {0,0},fact_bound[2];
90   PetscReal          rnorm,fnorm;
91 
92   PetscFunctionBegin;
93   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
94   ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr);
95   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
96   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
97   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
98   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
99 
100   rnorm = 0.0;
101   for (i=0; i<n; i++) {
102     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]);
103     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++;
104     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++;
105     else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Can never get here");
106   }
107 
108   for (i=0; i<n; i++) {
109     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++;
110     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++;
111   }
112   ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
113   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
114   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
115   ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
116   ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
117   ierr = MPI_Allreduce(act,fact,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
118   ierr = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
119   fnorm = PetscSqrtReal(fnorm);
120 
121   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
122   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);
123 
124   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
125   PetscFunctionReturn(0);
126 }
127 
128 /*
129      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
130     || 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
131     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
132     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
133 */
134 #undef __FUNCT__
135 #define __FUNCT__ "SNESVICheckLocalMin_Private"
136 PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
137 {
138   PetscReal      a1;
139   PetscErrorCode ierr;
140   PetscBool     hastranspose;
141 
142   PetscFunctionBegin;
143   *ismin = PETSC_FALSE;
144   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
145   if (hastranspose) {
146     /* Compute || J^T F|| */
147     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
148     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
149     ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
150     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
151   } else {
152     Vec         work;
153     PetscScalar result;
154     PetscReal   wnorm;
155 
156     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
157     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
158     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
159     ierr = MatMult(A,W,work);CHKERRQ(ierr);
160     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
161     ierr = VecDestroy(&work);CHKERRQ(ierr);
162     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
163     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
164     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
165   }
166   PetscFunctionReturn(0);
167 }
168 
169 /*
170      Checks if J^T(F - J*X) = 0
171 */
172 #undef __FUNCT__
173 #define __FUNCT__ "SNESVICheckResidual_Private"
174 PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
175 {
176   PetscReal      a1,a2;
177   PetscErrorCode ierr;
178   PetscBool     hastranspose;
179 
180   PetscFunctionBegin;
181   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
182   if (hastranspose) {
183     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
184     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
185 
186     /* Compute || J^T W|| */
187     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
188     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
189     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
190     if (a1 != 0.0) {
191       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
192     }
193   }
194   PetscFunctionReturn(0);
195 }
196 
197 /*
198   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
199 
200   Notes:
201   The convergence criterion currently implemented is
202   merit < abstol
203   merit < rtol*merit_initial
204 */
205 #undef __FUNCT__
206 #define __FUNCT__ "SNESDefaultConverged_VI"
207 PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
208 {
209   PetscErrorCode ierr;
210 
211   PetscFunctionBegin;
212   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
213   PetscValidPointer(reason,6);
214 
215   *reason = SNES_CONVERGED_ITERATING;
216 
217   if (!it) {
218     /* set parameter for default relative tolerance convergence test */
219     snes->ttol = fnorm*snes->rtol;
220   }
221   if (fnorm != fnorm) {
222     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
223     *reason = SNES_DIVERGED_FNORM_NAN;
224   } else if (fnorm < snes->abstol) {
225     ierr = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
226     *reason = SNES_CONVERGED_FNORM_ABS;
227   } else if (snes->nfuncs >= snes->max_funcs) {
228     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
229     *reason = SNES_DIVERGED_FUNCTION_COUNT;
230   }
231 
232   if (it && !*reason) {
233     if (fnorm < snes->ttol) {
234       ierr = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
235       *reason = SNES_CONVERGED_FNORM_RELATIVE;
236     }
237   }
238   PetscFunctionReturn(0);
239 }
240 
241 
242 /* -------------------------------------------------------------------------- */
243 /*
244    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
245 
246    Input Parameters:
247 .  SNES - nonlinear solver context
248 
249    Output Parameters:
250 .  X - Bound projected X
251 
252 */
253 
254 #undef __FUNCT__
255 #define __FUNCT__ "SNESVIProjectOntoBounds"
256 PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
257 {
258   PetscErrorCode    ierr;
259   const PetscScalar *xl,*xu;
260   PetscScalar       *x;
261   PetscInt          i,n;
262 
263   PetscFunctionBegin;
264   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
265   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
266   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
267   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
268 
269   for(i = 0;i<n;i++) {
270     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
271     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
272   }
273   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
274   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
275   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
276   PetscFunctionReturn(0);
277 }
278 
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "SNESVIGetActiveSetIS"
282 /*
283    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
284 
285    Input parameter
286 .  snes - the SNES context
287 .  X    - the snes solution vector
288 .  F    - the nonlinear function vector
289 
290    Output parameter
291 .  ISact - active set index set
292  */
293 PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS* ISact)
294 {
295   PetscErrorCode   ierr;
296   Vec               Xl=snes->xl,Xu=snes->xu;
297   const PetscScalar *x,*f,*xl,*xu;
298   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
299 
300   PetscFunctionBegin;
301   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
302   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
303   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
304   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
305   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
306   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
307   /* Compute active set size */
308   for (i=0; i < nlocal;i++) {
309     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++;
310   }
311 
312   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
313 
314   /* Set active set indices */
315   for(i=0; i < nlocal; i++) {
316     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;
317   }
318 
319    /* Create active set IS */
320   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
321 
322   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
323   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
324   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
325   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
326   PetscFunctionReturn(0);
327 }
328 
329 #undef __FUNCT__
330 #define __FUNCT__ "SNESVICreateIndexSets_RS"
331 PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS* ISact,IS* ISinact)
332 {
333   PetscErrorCode     ierr;
334 
335   PetscFunctionBegin;
336   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
337   ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr);
338   PetscFunctionReturn(0);
339 }
340 
341 #undef __FUNCT__
342 #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
343 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X,PetscReal *fnorm)
344 {
345   PetscErrorCode    ierr;
346   const PetscScalar *x,*xl,*xu,*f;
347   PetscInt          i,n;
348   PetscReal         rnorm;
349 
350   PetscFunctionBegin;
351   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
352   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
353   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
354   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
355   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
356   rnorm = 0.0;
357   for (i=0; i<n; i++) {
358     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]);
359   }
360   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
361   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
362   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
363   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
364   ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
365   *fnorm = PetscSqrtReal(*fnorm);
366   PetscFunctionReturn(0);
367 }
368 
369 
370 /* -------------------------------------------------------------------------- */
371 /*
372    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
373    of the SNESVI nonlinear solver.
374 
375    Input Parameter:
376 .  snes - the SNES context
377 .  x - the solution vector
378 
379    Application Interface Routine: SNESSetUp()
380 
381    Notes:
382    For basic use of the SNES solvers, the user need not explicitly call
383    SNESSetUp(), since these actions will automatically occur during
384    the call to SNESSolve().
385  */
386 #undef __FUNCT__
387 #define __FUNCT__ "SNESSetUp_VI"
388 PetscErrorCode SNESSetUp_VI(SNES snes)
389 {
390   PetscErrorCode ierr;
391   PetscInt       i_start[3],i_end[3];
392 
393   PetscFunctionBegin;
394 
395   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
396 
397   if (snes->ops->computevariablebounds) {
398     if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
399     if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
400     ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
401   } else if (!snes->xl && !snes->xu) {
402     /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
403     ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr);
404     ierr = VecSet(snes->xl,SNES_VI_NINF);CHKERRQ(ierr);
405     ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr);
406     ierr = VecSet(snes->xu,SNES_VI_INF);CHKERRQ(ierr);
407   } else {
408     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
409     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
410     ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr);
411     ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr);
412     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]))
413       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.");
414   }
415 
416   PetscFunctionReturn(0);
417 }
418 /* -------------------------------------------------------------------------- */
419 #undef __FUNCT__
420 #define __FUNCT__ "SNESReset_VI"
421 PetscErrorCode SNESReset_VI(SNES snes)
422 {
423   PetscErrorCode ierr;
424 
425   PetscFunctionBegin;
426   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
427   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
428   PetscFunctionReturn(0);
429 }
430 
431 /*
432    SNESDestroy_VI - Destroys the private SNES_VI context that was created
433    with SNESCreate_VI().
434 
435    Input Parameter:
436 .  snes - the SNES context
437 
438    Application Interface Routine: SNESDestroy()
439  */
440 #undef __FUNCT__
441 #define __FUNCT__ "SNESDestroy_VI"
442 PetscErrorCode SNESDestroy_VI(SNES snes)
443 {
444   PetscErrorCode ierr;
445 
446   PetscFunctionBegin;
447   ierr = PetscFree(snes->data);CHKERRQ(ierr);
448 
449   /* clear composed functions */
450   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
451   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
452   PetscFunctionReturn(0);
453 }
454 
455 /* -------------------------------------------------------------------------- */
456 
457 /*
458 
459       These line searches are common for all the VI solvers
460 */
461 extern PetscErrorCode SNESSolve_VISS(SNES);
462 
463 #undef __FUNCT__
464 #define __FUNCT__ "SNESLineSearchNo_VI"
465 
466 /*
467   This routine does not actually do a line search but it takes a full newton
468   step while ensuring that the new iterates remain within the constraints.
469 
470 */
471 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)
472 {
473   PetscErrorCode ierr;
474   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
475 
476   PetscFunctionBegin;
477   *flag = PETSC_TRUE;
478   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
479   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
480   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
481   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
482   if (snes->ops->postcheckstep) {
483    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
484   }
485   if (changed_y) {
486     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
487     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
488   }
489   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
490   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
491   if (!snes->domainerror) {
492     if (snes->ops->solve != SNESSolve_VISS) {
493        ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
494     } else {
495       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
496     }
497     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
498   }
499   if (snes->ls_monitor) {
500     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
501     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
502     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
503   }
504   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
505   PetscFunctionReturn(0);
506 }
507 
508 /* -------------------------------------------------------------------------- */
509 #undef __FUNCT__
510 #define __FUNCT__ "SNESLineSearchNoNorms_VI"
511 
512 /*
513   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
514 */
515 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)
516 {
517   PetscErrorCode ierr;
518   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
519 
520   PetscFunctionBegin;
521   *flag = PETSC_TRUE;
522   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
523   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
524   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
525   if (snes->ops->postcheckstep) {
526    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
527   }
528   if (changed_y) {
529     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
530     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
531   }
532 
533   /* don't evaluate function the last time through */
534   if (snes->iter < snes->max_its-1) {
535     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
536   }
537   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
538   PetscFunctionReturn(0);
539 }
540 
541 /* -------------------------------------------------------------------------- */
542 
543 #undef __FUNCT__
544 #define __FUNCT__ "SNESLineSearchCubic_VI"
545 /*
546   This routine implements a cubic line search while doing a projection on the variable bounds
547 */
548 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)
549 {
550   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
551   PetscReal      minlambda,lambda,lambdatemp;
552 #if defined(PETSC_USE_COMPLEX)
553   PetscScalar    cinitslope;
554 #endif
555   PetscErrorCode ierr;
556   PetscInt       count;
557   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
558   MPI_Comm       comm;
559 
560   PetscFunctionBegin;
561   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
562   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
563   *flag   = PETSC_TRUE;
564 
565   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
566   if (*ynorm == 0.0) {
567     if (snes->ls_monitor) {
568       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
569       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
570       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
571     }
572     *gnorm = fnorm;
573     ierr   = VecCopy(x,w);CHKERRQ(ierr);
574     ierr   = VecCopy(f,g);CHKERRQ(ierr);
575     *flag  = PETSC_FALSE;
576     goto theend1;
577   }
578   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
579     if (snes->ls_monitor) {
580       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
581       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Scaling step by %g old ynorm %g\n",(double)(snes->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr);
582       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
583     }
584     ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
585     *ynorm = snes->maxstep;
586   }
587   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
588   minlambda = snes->steptol/rellength;
589   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
590 #if defined(PETSC_USE_COMPLEX)
591   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
592   initslope = PetscRealPart(cinitslope);
593 #else
594   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
595 #endif
596   if (initslope > 0.0)  initslope = -initslope;
597   if (initslope == 0.0) initslope = -1.0;
598 
599   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
600   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
601   if (snes->nfuncs >= snes->max_funcs) {
602     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
603     *flag = PETSC_FALSE;
604     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
605     goto theend1;
606   }
607   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
608   if (snes->ops->solve != SNESSolve_VISS) {
609     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
610   } else {
611     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
612   }
613   if (snes->domainerror) {
614     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
615     PetscFunctionReturn(0);
616   }
617   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
618   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);
619   if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* Sufficient reduction */
620     if (snes->ls_monitor) {
621       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
622       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr);
623       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
624     }
625     goto theend1;
626   }
627 
628   /* Fit points with quadratic */
629   lambda     = 1.0;
630   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
631   lambdaprev = lambda;
632   gnormprev  = *gnorm;
633   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
634   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
635   else                         lambda = lambdatemp;
636 
637   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
638   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
639   if (snes->nfuncs >= snes->max_funcs) {
640     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
641     *flag = PETSC_FALSE;
642     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
643     goto theend1;
644   }
645   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
646   if (snes->ops->solve != SNESSolve_VISS) {
647     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
648   } else {
649     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
650   }
651   if (snes->domainerror) {
652     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
653     PetscFunctionReturn(0);
654   }
655   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
656   if (snes->ls_monitor) {
657     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
658     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: gnorm after quadratic fit %g\n",(double)(*gnorm));CHKERRQ(ierr);
659     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
660   }
661   if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* sufficient reduction */
662     if (snes->ls_monitor) {
663       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
664       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
665       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
666     }
667     goto theend1;
668   }
669 
670   /* Fit points with cubic */
671   count = 1;
672   while (PETSC_TRUE) {
673     if (lambda <= minlambda) {
674       if (snes->ls_monitor) {
675         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
676  	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
677 	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);
678         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
679       }
680       *flag = PETSC_FALSE;
681       break;
682     }
683     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
684     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
685     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
686     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
687     d  = b*b - 3*a*initslope;
688     if (d < 0.0) d = 0.0;
689     if (a == 0.0) {
690       lambdatemp = -initslope/(2.0*b);
691     } else {
692       lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
693     }
694     lambdaprev = lambda;
695     gnormprev  = *gnorm;
696     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
697     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
698     else                         lambda     = lambdatemp;
699     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
700     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
701     if (snes->nfuncs >= snes->max_funcs) {
702       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
703       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);
704       *flag = PETSC_FALSE;
705       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
706       break;
707     }
708     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
709     if (snes->ops->solve != SNESSolve_VISS) {
710       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
711     } else {
712       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
713     }
714     if (snes->domainerror) {
715       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
716       PetscFunctionReturn(0);
717     }
718     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
719     if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* is reduction enough? */
720       if (snes->ls_monitor) {
721 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %g lambda=%18.16e\n",(double)(*gnorm),(double)lambda);CHKERRQ(ierr);
722       }
723       break;
724     } else {
725       if (snes->ls_monitor) {
726         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %g lambda=%18.16e\n",(double)(*gnorm),(double)lambda);CHKERRQ(ierr);
727       }
728     }
729     count++;
730   }
731   theend1:
732   /* Optional user-defined check for line search step validity */
733   if (snes->ops->postcheckstep && *flag) {
734     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
735     if (changed_y) {
736       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
737       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
738     }
739     if (changed_y || changed_w) { /* recompute the function if the step has changed */
740       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
741       if (snes->ops->solve != SNESSolve_VISS) {
742         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
743       } else {
744         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
745       }
746       if (snes->domainerror) {
747         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
748         PetscFunctionReturn(0);
749       }
750       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
751       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
752       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
753     }
754   }
755   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
756   PetscFunctionReturn(0);
757 }
758 
759 #undef __FUNCT__
760 #define __FUNCT__ "SNESLineSearchQuadratic_VI"
761 /*
762   This routine does a quadratic line search while keeping the iterates within the variable bounds
763 */
764 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)
765 {
766   /*
767      Note that for line search purposes we work with with the related
768      minimization problem:
769         min  z(x):  R^n -> R,
770      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
771    */
772   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
773 #if defined(PETSC_USE_COMPLEX)
774   PetscScalar    cinitslope;
775 #endif
776   PetscErrorCode ierr;
777   PetscInt       count;
778   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
779 
780   PetscFunctionBegin;
781   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
782   *flag   = PETSC_TRUE;
783 
784   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
785   if (*ynorm == 0.0) {
786     if (snes->ls_monitor) {
787       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
788       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
789       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
790     }
791     *gnorm = fnorm;
792     ierr   = VecCopy(x,w);CHKERRQ(ierr);
793     ierr   = VecCopy(f,g);CHKERRQ(ierr);
794     *flag  = PETSC_FALSE;
795     goto theend2;
796   }
797   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
798     ierr   = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
799     *ynorm = snes->maxstep;
800   }
801   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
802   minlambda = snes->steptol/rellength;
803   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
804 #if defined(PETSC_USE_COMPLEX)
805   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
806   initslope = PetscRealPart(cinitslope);
807 #else
808   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
809 #endif
810   if (initslope > 0.0)  initslope = -initslope;
811   if (initslope == 0.0) initslope = -1.0;
812   ierr = PetscInfo1(snes,"Initslope %g \n",(double)initslope);CHKERRQ(ierr);
813 
814   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
815   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
816   if (snes->nfuncs >= snes->max_funcs) {
817     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
818     *flag = PETSC_FALSE;
819     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
820     goto theend2;
821   }
822   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
823   if (snes->ops->solve != SNESSolve_VISS) {
824     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
825   } else {
826     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
827   }
828   if (snes->domainerror) {
829     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
830     PetscFunctionReturn(0);
831   }
832   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
833   if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* Sufficient reduction */
834     if (snes->ls_monitor) {
835       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
836       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr);
837       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
838     }
839     goto theend2;
840   }
841 
842   /* Fit points with quadratic */
843   lambda = 1.0;
844   count = 1;
845   while (PETSC_TRUE) {
846     if (lambda <= minlambda) { /* bad luck; use full step */
847       if (snes->ls_monitor) {
848         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
849         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
850         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);
851         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
852       }
853       ierr = VecCopy(x,w);CHKERRQ(ierr);
854       *flag = PETSC_FALSE;
855       break;
856     }
857     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
858     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
859     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
860     else                         lambda     = lambdatemp;
861 
862     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
863     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
864     if (snes->nfuncs >= snes->max_funcs) {
865       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
866       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);
867       *flag = PETSC_FALSE;
868       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
869       break;
870     }
871     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
872     if (snes->domainerror) {
873       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
874       PetscFunctionReturn(0);
875     }
876     if (snes->ops->solve != SNESSolve_VISS) {
877       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
878     } else {
879       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
880     }
881     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
882     if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* sufficient reduction */
883       if (snes->ls_monitor) {
884         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
885         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Quadratically determined step, lambda=%g\n",(double)lambda);CHKERRQ(ierr);
886         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
887       }
888       break;
889     }
890     count++;
891   }
892   theend2:
893   /* Optional user-defined check for line search step validity */
894   if (snes->ops->postcheckstep) {
895     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
896     if (changed_y) {
897       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
898       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
899     }
900     if (changed_y || changed_w) { /* recompute the function if the step has changed */
901       ierr = SNESComputeFunction(snes,w,g);
902       if (snes->domainerror) {
903         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
904         PetscFunctionReturn(0);
905       }
906       if (snes->ops->solve != SNESSolve_VISS) {
907         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
908       } else {
909         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
910       }
911 
912       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
913       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
914       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
915     }
916   }
917   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
918   PetscFunctionReturn(0);
919 }
920 
921 
922 #undef __FUNCT__
923 #define __FUNCT__ "SNESVISetVariableBounds"
924 /*@
925    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
926 
927    Input Parameters:
928 .  snes - the SNES context.
929 .  xl   - lower bound.
930 .  xu   - upper bound.
931 
932    Notes:
933    If this routine is not called then the lower and upper bounds are set to
934    SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp().
935 
936    Level: advanced
937 
938 @*/
939 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
940 {
941   PetscErrorCode    ierr;
942   const PetscScalar *xxl,*xxu;
943   PetscInt          i,n, cnt = 0;
944 
945   PetscFunctionBegin;
946   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
947   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
948   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
949   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
950   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);
951   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);
952   ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);
953   ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
954   ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
955   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
956   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
957   snes->xl = xl;
958   snes->xu = xu;
959   ierr = VecGetLocalSize(xl,&n);CHKERRQ(ierr);
960   ierr = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr);
961   ierr = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr);
962   for (i=0; i<n; i++) {
963     cnt += ((xxl[i] != SNES_VI_NINF) || (xxu[i] != SNES_VI_INF));
964   }
965   ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
966   ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr);
967   ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr);
968   PetscFunctionReturn(0);
969 }
970 
971 
972 EXTERN_C_BEGIN
973 #undef __FUNCT__
974 #define __FUNCT__ "SNESLineSearchSetType_VI"
975 PetscErrorCode  SNESLineSearchSetType_VI(SNES snes, SNESLineSearchType type)
976 {
977   PetscErrorCode ierr;
978   PetscFunctionBegin;
979 
980   switch (type) {
981   case SNES_LS_BASIC:
982     ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
983     break;
984   case SNES_LS_BASIC_NONORMS:
985     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
986     break;
987   case SNES_LS_QUADRATIC:
988     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
989     break;
990   case SNES_LS_CUBIC:
991     ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
992     break;
993   default:
994     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
995     break;
996   }
997   snes->ls_type = type;
998   PetscFunctionReturn(0);
999 }
1000 EXTERN_C_END
1001 
1002 #undef __FUNCT__
1003 #define __FUNCT__ "SNESSetFromOptions_VI"
1004 PetscErrorCode SNESSetFromOptions_VI(SNES snes)
1005 {
1006   PetscErrorCode  ierr;
1007   PetscBool       flg;
1008 
1009   PetscFunctionBegin;
1010   ierr = PetscOptionsHead("SNES VI options");CHKERRQ(ierr);
1011   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
1012   if (flg) {
1013     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
1014   }
1015   ierr = PetscOptionsTail();CHKERRQ(ierr);
1016   PetscFunctionReturn(0);
1017 }
1018