xref: /petsc/src/snes/impls/vi/vi.c (revision abcba3417629b6ec761a96e9f095a7df6d3bf4ed)
1 #define PETSCSNES_DLL
2 
3 #include "../src/snes/impls/vi/viimpl.h" /*I "petscsnes.h" I*/
4 #include "../include/private/kspimpl.h"
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "SNESMonitorVI"
8 PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
9 {
10   PetscErrorCode          ierr;
11   SNES_VI                 *vi = (SNES_VI*)snes->data;
12   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy;
13   const PetscScalar       *x,*xl,*xu,*f;
14   PetscInt                i,n,act = 0;
15   PetscReal               rnorm,fnorm;
16 
17   PetscFunctionBegin;
18   if (!dummy) {
19     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
20   }
21   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
22   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
23   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
24   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
25   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
26 
27   rnorm = 0.0;
28   for (i=0; i<n; i++) {
29     if (((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) rnorm += f[i]*f[i];
30     else act++;
31   }
32   ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
33   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
34   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
35   ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
36   ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
37   fnorm = sqrt(fnorm);
38   ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES VI Function norm %14.12e Active constraints %D\n",its,fnorm,act);CHKERRQ(ierr);
39   if (!dummy) {
40     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
41   }
42   PetscFunctionReturn(0);
43 }
44 
45 /*
46      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
47     || 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
48     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
49     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
50 */
51 #undef __FUNCT__
52 #define __FUNCT__ "SNESVICheckLocalMin_Private"
53 PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
54 {
55   PetscReal      a1;
56   PetscErrorCode ierr;
57   PetscBool     hastranspose;
58 
59   PetscFunctionBegin;
60   *ismin = PETSC_FALSE;
61   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
62   if (hastranspose) {
63     /* Compute || J^T F|| */
64     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
65     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
66     ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr);
67     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
68   } else {
69     Vec         work;
70     PetscScalar result;
71     PetscReal   wnorm;
72 
73     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
74     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
75     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
76     ierr = MatMult(A,W,work);CHKERRQ(ierr);
77     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
78     ierr = VecDestroy(work);CHKERRQ(ierr);
79     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
80     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr);
81     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
82   }
83   PetscFunctionReturn(0);
84 }
85 
86 /*
87      Checks if J^T(F - J*X) = 0
88 */
89 #undef __FUNCT__
90 #define __FUNCT__ "SNESVICheckResidual_Private"
91 PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
92 {
93   PetscReal      a1,a2;
94   PetscErrorCode ierr;
95   PetscBool     hastranspose;
96 
97   PetscFunctionBegin;
98   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
99   if (hastranspose) {
100     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
101     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
102 
103     /* Compute || J^T W|| */
104     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
105     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
106     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
107     if (a1 != 0.0) {
108       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr);
109     }
110   }
111   PetscFunctionReturn(0);
112 }
113 
114 /*
115   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
116 
117   Notes:
118   The convergence criterion currently implemented is
119   merit < abstol
120   merit < rtol*merit_initial
121 */
122 #undef __FUNCT__
123 #define __FUNCT__ "SNESDefaultConverged_VI"
124 PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
125 {
126   PetscErrorCode ierr;
127 
128   PetscFunctionBegin;
129   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
130   PetscValidPointer(reason,6);
131 
132   *reason = SNES_CONVERGED_ITERATING;
133 
134   if (!it) {
135     /* set parameter for default relative tolerance convergence test */
136     snes->ttol = fnorm*snes->rtol;
137   }
138   if (fnorm != fnorm) {
139     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
140     *reason = SNES_DIVERGED_FNORM_NAN;
141   } else if (fnorm < snes->abstol) {
142     ierr = PetscInfo2(snes,"Converged due to function norm %G < %G\n",fnorm,snes->abstol);CHKERRQ(ierr);
143     *reason = SNES_CONVERGED_FNORM_ABS;
144   } else if (snes->nfuncs >= snes->max_funcs) {
145     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
146     *reason = SNES_DIVERGED_FUNCTION_COUNT;
147   }
148 
149   if (it && !*reason) {
150     if (fnorm < snes->ttol) {
151       ierr = PetscInfo2(snes,"Converged due to function norm %G < %G (relative tolerance)\n",fnorm,snes->ttol);CHKERRQ(ierr);
152       *reason = SNES_CONVERGED_FNORM_RELATIVE;
153     }
154   }
155   PetscFunctionReturn(0);
156 }
157 
158 /*
159   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
160 
161   Input Parameter:
162 . phi - the semismooth function
163 
164   Output Parameter:
165 . merit - the merit function
166 . phinorm - ||phi||
167 
168   Notes:
169   The merit function for the mixed complementarity problem is defined as
170      merit = 0.5*phi^T*phi
171 */
172 #undef __FUNCT__
173 #define __FUNCT__ "SNESVIComputeMeritFunction"
174 static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm)
175 {
176   PetscErrorCode ierr;
177 
178   PetscFunctionBegin;
179   ierr = VecNormBegin(phi,NORM_2,phinorm);
180   ierr = VecNormEnd(phi,NORM_2,phinorm);
181 
182   *merit = 0.5*(*phinorm)*(*phinorm);
183   PetscFunctionReturn(0);
184 }
185 
186 PETSC_STATIC_INLINE PetscScalar Phi(PetscScalar a,PetscScalar b)
187 {
188   return a + b - sqrt(a*a + b*b);
189 }
190 
191 PETSC_STATIC_INLINE PetscScalar DPhi(PetscScalar a,PetscScalar b)
192 {
193   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return  1.0 - a/ sqrt(a*a + b*b);
194   else return .5;
195 }
196 
197 /*
198    SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form.
199 
200    Input Parameters:
201 .  snes - the SNES context
202 .  x - current iterate
203 .  functx - user defined function context
204 
205    Output Parameters:
206 .  phi - Semismooth function
207 
208 */
209 #undef __FUNCT__
210 #define __FUNCT__ "SNESVIComputeFunction"
211 static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx)
212 {
213   PetscErrorCode  ierr;
214   SNES_VI         *vi = (SNES_VI*)snes->data;
215   Vec             Xl = vi->xl,Xu = vi->xu,F = snes->vec_func;
216   PetscScalar     *phi_arr,*x_arr,*f_arr,*l,*u;
217   PetscInt        i,nlocal;
218 
219   PetscFunctionBegin;
220   ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr);
221   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
222   ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr);
223   ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr);
224   ierr = VecGetArray(Xl,&l);CHKERRQ(ierr);
225   ierr = VecGetArray(Xu,&u);CHKERRQ(ierr);
226   ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr);
227 
228   for (i=0;i < nlocal;i++) {
229     if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { /* no constraints on variable */
230       phi_arr[i] = f_arr[i];
231     } else if (l[i] <= PETSC_VI_NINF) {                      /* upper bound on variable only */
232       phi_arr[i] = -Phi(u[i] - x_arr[i],-f_arr[i]);
233     } else if (u[i] >= PETSC_VI_INF) {                       /* lower bound on variable only */
234       phi_arr[i] = Phi(x_arr[i] - l[i],f_arr[i]);
235     } else if (l[i] == u[i]) {
236       phi_arr[i] = l[i] - x_arr[i];
237     } else {                                                /* both bounds on variable */
238       phi_arr[i] = Phi(x_arr[i] - l[i],-Phi(u[i] - x_arr[i],-f_arr[i]));
239     }
240   }
241 
242   ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr);
243   ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr);
244   ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr);
245   ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr);
246   ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr);
247   PetscFunctionReturn(0);
248 }
249 
250 /*
251    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
252                                           the semismooth jacobian.
253 */
254 #undef __FUNCT__
255 #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors"
256 PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db)
257 {
258   PetscErrorCode ierr;
259   SNES_VI      *vi = (SNES_VI*)snes->data;
260   PetscScalar    *l,*u,*x,*f,*da,*db,da1,da2,db1,db2;
261   PetscInt       i,nlocal;
262 
263   PetscFunctionBegin;
264 
265   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
266   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
267   ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr);
268   ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr);
269   ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
270   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
271   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
272 
273   for (i=0;i< nlocal;i++) {
274     if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) {/* no constraints on variable */
275       da[i] = 0;
276       db[i] = 1;
277     } else if (l[i] <= PETSC_VI_NINF) {                     /* upper bound on variable only */
278       da[i] = DPhi(u[i] - x[i], -f[i]);
279       db[i] = DPhi(-f[i],u[i] - x[i]);
280     } else if (u[i] >= PETSC_VI_INF) {                      /* lower bound on variable only */
281       da[i] = DPhi(x[i] - l[i], f[i]);
282       db[i] = DPhi(f[i],x[i] - l[i]);
283     } else if (l[i] == u[i]) {                              /* fixed variable */
284       da[i] = 1;
285       db[i] = 0;
286     } else {                                                /* upper and lower bounds on variable */
287       da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
288       db1 = DPhi(-Phi(u[i] - x[i], -f[i]),x[i] - l[i]);
289       da2 = DPhi(u[i] - x[i], -f[i]);
290       db2 = DPhi(-f[i],u[i] - x[i]);
291       da[i] = da1 + db1*da2;
292       db[i] = db1*db2;
293     }
294   }
295 
296   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
297   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
298   ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr);
299   ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr);
300   ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
301   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
302   PetscFunctionReturn(0);
303 }
304 
305 /*
306    SNESVIComputeJacobian - Computes the jacobian of the semismooth function.The Jacobian for the semismooth function is an element of the B-subdifferential of the Fischer-Burmeister function for complementarity problems.
307 
308    Input Parameters:
309 .  Da       - Diagonal shift vector for the semismooth jacobian.
310 .  Db       - Row scaling vector for the semismooth jacobian.
311 
312    Output Parameters:
313 .  jac      - semismooth jacobian
314 .  jac_pre  - optional preconditioning matrix
315 
316    Notes:
317    The semismooth jacobian matrix is given by
318    jac = Da + Db*jacfun
319    where Db is the row scaling matrix stored as a vector,
320          Da is the diagonal perturbation matrix stored as a vector
321    and   jacfun is the jacobian of the original nonlinear function.
322 */
323 #undef __FUNCT__
324 #define __FUNCT__ "SNESVIComputeJacobian"
325 PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db)
326 {
327   PetscErrorCode ierr;
328 
329   /* Do row scaling  and add diagonal perturbation */
330   ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr);
331   ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr);
332   if (jac != jac_pre) { /* If jac and jac_pre are different */
333     ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL);
334     ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr);
335   }
336   PetscFunctionReturn(0);
337 }
338 
339 /*
340    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
341 
342    Input Parameters:
343    phi - semismooth function.
344    H   - semismooth jacobian
345 
346    Output Parameters:
347    dpsi - merit function gradient
348 
349    Notes:
350   The merit function gradient is computed as follows
351         dpsi = H^T*phi
352 */
353 #undef __FUNCT__
354 #define __FUNCT__ "SNESVIComputeMeritFunctionGradient"
355 PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi)
356 {
357   PetscErrorCode ierr;
358 
359   PetscFunctionBegin;
360   ierr = MatMultTranspose(H,phi,dpsi);
361   PetscFunctionReturn(0);
362 }
363 
364 /* -------------------------------------------------------------------------- */
365 /*
366    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
367 
368    Input Parameters:
369 .  SNES - nonlinear solver context
370 
371    Output Parameters:
372 .  X - Bound projected X
373 
374 */
375 
376 #undef __FUNCT__
377 #define __FUNCT__ "SNESVIProjectOntoBounds"
378 PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
379 {
380   PetscErrorCode    ierr;
381   SNES_VI           *vi = (SNES_VI*)snes->data;
382   const PetscScalar *xl,*xu;
383   PetscScalar       *x;
384   PetscInt          i,n;
385 
386   PetscFunctionBegin;
387   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
388   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
389   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
390   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
391 
392   for(i = 0;i<n;i++) {
393     if (x[i] < xl[i]) x[i] = xl[i];
394     else if (x[i] > xu[i]) x[i] = xu[i];
395   }
396   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
397   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
398   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
399   PetscFunctionReturn(0);
400 }
401 
402 /*  --------------------------------------------------------------------
403 
404      This file implements a semismooth truncated Newton method with a line search,
405      for solving a system of nonlinear equations in complementarity form, using the KSP, Vec,
406      and Mat interfaces for linear solvers, vectors, and matrices,
407      respectively.
408 
409      The following basic routines are required for each nonlinear solver:
410           SNESCreate_XXX()          - Creates a nonlinear solver context
411           SNESSetFromOptions_XXX()  - Sets runtime options
412           SNESSolve_XXX()           - Solves the nonlinear system
413           SNESDestroy_XXX()         - Destroys the nonlinear solver context
414      The suffix "_XXX" denotes a particular implementation, in this case
415      we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving
416      systems of nonlinear equations with a line search (LS) method.
417      These routines are actually called via the common user interface
418      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
419      SNESDestroy(), so the application code interface remains identical
420      for all nonlinear solvers.
421 
422      Another key routine is:
423           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
424      by setting data structures and options.   The interface routine SNESSetUp()
425      is not usually called directly by the user, but instead is called by
426      SNESSolve() if necessary.
427 
428      Additional basic routines are:
429           SNESView_XXX()            - Prints details of runtime options that
430                                       have actually been used.
431      These are called by application codes via the interface routines
432      SNESView().
433 
434      The various types of solvers (preconditioners, Krylov subspace methods,
435      nonlinear solvers, timesteppers) are all organized similarly, so the
436      above description applies to these categories also.
437 
438     -------------------------------------------------------------------- */
439 /*
440    SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton
441    method using a line search.
442 
443    Input Parameters:
444 .  snes - the SNES context
445 
446    Output Parameter:
447 .  outits - number of iterations until termination
448 
449    Application Interface Routine: SNESSolve()
450 
451    Notes:
452    This implements essentially a semismooth Newton method with a
453    line search. The default line search does not do any line seach
454    but rather takes a full newton step.
455 */
456 #undef __FUNCT__
457 #define __FUNCT__ "SNESSolveVI_SS"
458 PetscErrorCode SNESSolveVI_SS(SNES snes)
459 {
460   SNES_VI            *vi = (SNES_VI*)snes->data;
461   PetscErrorCode     ierr;
462   PetscInt           maxits,i,lits;
463   PetscBool          lssucceed;
464   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
465   PetscReal          gnorm,xnorm=0,ynorm;
466   Vec                Y,X,F,G,W;
467   KSPConvergedReason kspreason;
468 
469   PetscFunctionBegin;
470   vi->computeuserfunction    = snes->ops->computefunction;
471   snes->ops->computefunction = SNESVIComputeFunction;
472 
473   snes->numFailures            = 0;
474   snes->numLinearSolveFailures = 0;
475   snes->reason                 = SNES_CONVERGED_ITERATING;
476 
477   maxits	= snes->max_its;	/* maximum number of iterations */
478   X		= snes->vec_sol;	/* solution vector */
479   F		= snes->vec_func;	/* residual vector */
480   Y		= snes->work[0];	/* work vectors */
481   G		= snes->work[1];
482   W		= snes->work[2];
483 
484   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
485   snes->iter = 0;
486   snes->norm = 0.0;
487   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
488 
489   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
490   ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr);
491   if (snes->domainerror) {
492     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
493     snes->ops->computefunction = vi->computeuserfunction;
494     PetscFunctionReturn(0);
495   }
496    /* Compute Merit function */
497   ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr);
498 
499   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
500   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
501   if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
502 
503   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
504   snes->norm = vi->phinorm;
505   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
506   SNESLogConvHistory(snes,vi->phinorm,0);
507   SNESMonitor(snes,0,vi->phinorm);
508 
509   /* set parameter for default relative tolerance convergence test */
510   snes->ttol = vi->phinorm*snes->rtol;
511   /* test convergence */
512   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
513   if (snes->reason) {
514     snes->ops->computefunction = vi->computeuserfunction;
515     PetscFunctionReturn(0);
516   }
517 
518   for (i=0; i<maxits; i++) {
519 
520     /* Call general purpose update function */
521     if (snes->ops->update) {
522       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
523     }
524 
525     /* Solve J Y = Phi, where J is the semismooth jacobian */
526     /* Get the nonlinear function jacobian */
527     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
528     /* Get the diagonal shift and row scaling vectors */
529     ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr);
530     /* Compute the semismooth jacobian */
531     ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr);
532     /* Compute the merit function gradient */
533     ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr);
534     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
535     ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr);
536     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
537 
538     if (kspreason < 0) {
539       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
540         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
541         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
542         break;
543       }
544     }
545     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
546     snes->linear_its += lits;
547     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
548     /*
549     if (vi->precheckstep) {
550       PetscBool changed_y = PETSC_FALSE;
551       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
552     }
553 
554     if (PetscLogPrintInfo){
555       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
556     }
557     */
558     /* Compute a (scaled) negative update in the line search routine:
559          Y <- X - lambda*Y
560        and evaluate G = function(Y) (depends on the line search).
561     */
562     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
563     ynorm = 1; gnorm = vi->phinorm;
564     ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
565     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
566     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
567     if (snes->domainerror) {
568       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
569       snes->ops->computefunction = vi->computeuserfunction;
570       PetscFunctionReturn(0);
571     }
572     if (!lssucceed) {
573       if (++snes->numFailures >= snes->maxFailures) {
574 	PetscBool ismin;
575         snes->reason = SNES_DIVERGED_LINE_SEARCH;
576         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
577         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
578         break;
579       }
580     }
581     /* Update function and solution vectors */
582     vi->phinorm = gnorm;
583     vi->merit = 0.5*vi->phinorm*vi->phinorm;
584     ierr = VecCopy(G,vi->phi);CHKERRQ(ierr);
585     ierr = VecCopy(W,X);CHKERRQ(ierr);
586     /* Monitor convergence */
587     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
588     snes->iter = i+1;
589     snes->norm = vi->phinorm;
590     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
591     SNESLogConvHistory(snes,snes->norm,lits);
592     SNESMonitor(snes,snes->iter,snes->norm);
593     /* Test for convergence, xnorm = || X || */
594     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
595     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
596     if (snes->reason) break;
597   }
598   if (i == maxits) {
599     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
600     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
601   }
602   snes->ops->computefunction = vi->computeuserfunction;
603   PetscFunctionReturn(0);
604 }
605 
606 #undef __FUNCT__
607 #define __FUNCT__ "SNESVIGetActiveSetIS"
608 /*
609    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
610 
611    Input parameter
612 .  snes - the SNES context
613 .  X    - the snes solution vector
614 .  F    - the nonlinear function vector
615 
616    Output parameter
617 .  ISact - active set index set
618  */
619 PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS* ISact)
620 {
621   PetscErrorCode   ierr;
622   SNES_VI          *vi = (SNES_VI*)snes->data;
623   Vec               Xl=vi->xl,Xu=vi->xu;
624   const PetscScalar *x,*f,*xl,*xu;
625   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
626 
627   PetscFunctionBegin;
628   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
629   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
630   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
631   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
632   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
633   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
634   /* Compute active set size */
635   for (i=0; i < nlocal;i++) {
636     if (!((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) nloc_isact++;
637   }
638 
639   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
640 
641   /* Set active set indices */
642   for(i=0; i < nlocal; i++) {
643     if (!((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) idx_act[i1++] = ilow+i;
644   }
645 
646    /* Create active set IS */
647   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
648 
649   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
650   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
651   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
652   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
653   PetscFunctionReturn(0);
654 }
655 
656 #undef __FUNCT__
657 #define __FUNCT__ "SNESVICreateIndexSets_RS"
658 PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS* ISact,IS* ISinact)
659 {
660   PetscErrorCode     ierr;
661 
662   PetscFunctionBegin;
663   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
664   ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr);
665   PetscFunctionReturn(0);
666 }
667 
668 /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */
669 #undef __FUNCT__
670 #define __FUNCT__ "SNESVICreateSubVectors"
671 PetscErrorCode SNESVICreateSubVectors(SNES snes,PetscInt n,Vec* newv)
672 {
673   PetscErrorCode ierr;
674   Vec            v;
675 
676   PetscFunctionBegin;
677   ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr);
678   ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr);
679   ierr = VecSetFromOptions(v);CHKERRQ(ierr);
680   *newv = v;
681 
682   PetscFunctionReturn(0);
683 }
684 
685 /* Resets the snes PC and KSP when the active set sizes change */
686 #undef __FUNCT__
687 #define __FUNCT__ "SNESVIResetPCandKSP"
688 PetscErrorCode SNESVIResetPCandKSP(SNES snes,Mat Amat,Mat Pmat)
689 {
690   PetscErrorCode         ierr;
691   KSP                    kspnew,snesksp;
692   PC                     pcnew;
693   const MatSolverPackage stype;
694 
695   PetscFunctionBegin;
696   /* The active and inactive set sizes have changed so need to create a new snes->ksp object */
697   ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr);
698   ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr);
699   /* Copy over snes->ksp info */
700   kspnew->pc_side = snesksp->pc_side;
701   kspnew->rtol    = snesksp->rtol;
702   kspnew->abstol    = snesksp->abstol;
703   kspnew->max_it  = snesksp->max_it;
704   ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr);
705   ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr);
706   ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr);
707   ierr = PCSetOperators(kspnew->pc,Amat,Pmat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
708   ierr = PCFactorGetMatSolverPackage(snesksp->pc,&stype);CHKERRQ(ierr);
709   ierr = PCFactorSetMatSolverPackage(kspnew->pc,stype);CHKERRQ(ierr);
710   ierr = KSPDestroy(snesksp);CHKERRQ(ierr);
711   snes->ksp = kspnew;
712   ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr);
713   ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr);
714   PetscFunctionReturn(0);
715 }
716 
717 
718 #undef __FUNCT__
719 #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
720 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X,PetscScalar *fnorm)
721 {
722   PetscErrorCode    ierr;
723   SNES_VI           *vi = (SNES_VI*)snes->data;
724   const PetscScalar *x,*xl,*xu,*f;
725   PetscInt          i,n;
726   PetscReal         rnorm;
727 
728   PetscFunctionBegin;
729   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
730   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
731   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
732   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
733   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
734   rnorm = 0.0;
735   for (i=0; i<n; i++) {
736     if (((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) rnorm += f[i]*f[i];
737   }
738   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
739   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
740   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
741   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
742   ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
743   *fnorm = sqrt(*fnorm);
744   PetscFunctionReturn(0);
745 }
746 
747 /* Variational Inequality solver using reduce space method. No semismooth algorithm is
748    implemented in this algorithm. It basically identifies the active variables and does
749    a linear solve on the inactive variables. */
750 #undef __FUNCT__
751 #define __FUNCT__ "SNESSolveVI_RS"
752 PetscErrorCode SNESSolveVI_RS(SNES snes)
753 {
754   SNES_VI          *vi = (SNES_VI*)snes->data;
755   PetscErrorCode    ierr;
756   PetscInt          maxits,i,lits;
757   PetscBool         lssucceed;
758   MatStructure      flg = DIFFERENT_NONZERO_PATTERN;
759   PetscReal         fnorm,gnorm,xnorm=0,ynorm;
760   Vec                Y,X,F,G,W;
761   KSPConvergedReason kspreason;
762 
763   PetscFunctionBegin;
764   snes->numFailures            = 0;
765   snes->numLinearSolveFailures = 0;
766   snes->reason                 = SNES_CONVERGED_ITERATING;
767 
768   maxits	= snes->max_its;	/* maximum number of iterations */
769   X		= snes->vec_sol;	/* solution vector */
770   F		= snes->vec_func;	/* residual vector */
771   Y		= snes->work[0];	/* work vectors */
772   G		= snes->work[1];
773   W		= snes->work[2];
774 
775   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
776   snes->iter = 0;
777   snes->norm = 0.0;
778   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
779 
780   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
781   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
782   if (snes->domainerror) {
783     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
784     PetscFunctionReturn(0);
785   }
786   ierr = SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);CHKERRQ(ierr);
787   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
788   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
789   if PetscIsInfOrNanReal(fnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
790 
791   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
792   snes->norm = fnorm;
793   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
794   SNESLogConvHistory(snes,fnorm,0);
795   SNESMonitor(snes,0,fnorm);
796 
797   /* set parameter for default relative tolerance convergence test */
798   snes->ttol = fnorm*snes->rtol;
799   /* test convergence */
800   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
801   if (snes->reason) PetscFunctionReturn(0);
802 
803   for (i=0; i<maxits; i++) {
804 
805     IS         IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
806     VecScatter scat_act,scat_inact;
807     PetscInt   nis_act,nis_inact;
808     Vec        Y_act,Y_inact,F_inact;
809     Mat        jac_inact_inact,prejac_inact_inact;
810     IS         keptrows;
811     PetscBool  isequal;
812 
813     /* Call general purpose update function */
814     if (snes->ops->update) {
815       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
816     }
817     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
818 
819     /* Create active and inactive index sets */
820     ierr = SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&IS_inact);CHKERRQ(ierr);
821 
822     /* Create inactive set submatrix */
823     ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
824     ierr = MatFindNonzeroRows(jac_inact_inact,&keptrows);CHKERRQ(ierr);
825     if (keptrows) {
826       PetscInt       cnt,*nrows,k;
827       const PetscInt *krows,*inact;
828 
829       ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr);
830       ierr = ISDestroy(IS_act);CHKERRQ(ierr);
831 
832       ierr = ISGetLocalSize(keptrows,&cnt);CHKERRQ(ierr);
833       ierr = ISGetIndices(keptrows,&krows);CHKERRQ(ierr);
834       ierr = ISGetIndices(IS_inact,&inact);CHKERRQ(ierr);
835       ierr = PetscMalloc(cnt*sizeof(PetscInt),&nrows);CHKERRQ(ierr);
836       for (k=0; k<cnt; k++) {
837         nrows[k] = inact[krows[k]];
838       }
839       ierr = ISRestoreIndices(keptrows,&krows);CHKERRQ(ierr);
840       ierr = ISRestoreIndices(IS_inact,&inact);CHKERRQ(ierr);
841       ierr = ISDestroy(keptrows);CHKERRQ(ierr);
842       ierr = ISDestroy(IS_inact);CHKERRQ(ierr);
843 
844       ierr = ISCreateGeneral(PETSC_COMM_SELF,cnt,nrows,PETSC_OWN_POINTER,&IS_inact);CHKERRQ(ierr);
845       ierr = ISComplement(IS_inact,0,F->map->n,&IS_act);CHKERRQ(ierr);
846       ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
847     }
848 
849     /* Get sizes of active and inactive sets */
850     ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr);
851     ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr);
852 
853     /* Create active and inactive set vectors */
854     ierr = SNESVICreateSubVectors(snes,nis_inact,&F_inact);CHKERRQ(ierr);
855     ierr = SNESVICreateSubVectors(snes,nis_act,&Y_act);CHKERRQ(ierr);
856     ierr = SNESVICreateSubVectors(snes,nis_inact,&Y_inact);CHKERRQ(ierr);
857 
858     /* Create scatter contexts */
859     ierr = VecScatterCreate(Y,IS_act,Y_act,PETSC_NULL,&scat_act);CHKERRQ(ierr);
860     ierr = VecScatterCreate(Y,IS_inact,Y_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr);
861 
862     /* Do a vec scatter to active and inactive set vectors */
863     ierr = VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
864     ierr = VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
865 
866     ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
867     ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
868 
869     ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
870     ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
871 
872     /* Active set direction = 0 */
873     ierr = VecSet(Y_act,0);CHKERRQ(ierr);
874     if (snes->jacobian != snes->jacobian_pre) {
875       ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr);
876     } else prejac_inact_inact = jac_inact_inact;
877 
878     ierr = ISEqual(vi->IS_inact_prev,IS_inact,&isequal);CHKERRQ(ierr);
879     if (!isequal) {
880       ierr = SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);CHKERRQ(ierr);
881     }
882 
883     ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr);
884     ierr = SNES_KSPSolve(snes,snes->ksp,F_inact,Y_inact);CHKERRQ(ierr);
885     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
886     if (kspreason < 0) {
887       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
888         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
889         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
890         break;
891       }
892      }
893 
894     ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
895     ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
896     ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
897     ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
898 
899     ierr = VecDestroy(F_inact);CHKERRQ(ierr);
900     ierr = VecDestroy(Y_act);CHKERRQ(ierr);
901     ierr = VecDestroy(Y_inact);CHKERRQ(ierr);
902     ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr);
903     ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr);
904     ierr = ISDestroy(IS_act);CHKERRQ(ierr);
905     if (!isequal) {
906       ierr = ISDestroy(vi->IS_inact_prev);CHKERRQ(ierr);
907       ierr = ISDuplicate(IS_inact,&vi->IS_inact_prev);CHKERRQ(ierr);
908     }
909     ierr = ISDestroy(IS_inact);CHKERRQ(ierr);
910     ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr);
911     if (snes->jacobian != snes->jacobian_pre) {
912       ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr);
913     }
914     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
915     snes->linear_its += lits;
916     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
917     /*
918     if (vi->precheckstep) {
919       PetscBool changed_y = PETSC_FALSE;
920       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
921     }
922 
923     if (PetscLogPrintInfo){
924       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
925     }
926     */
927     /* Compute a (scaled) negative update in the line search routine:
928          Y <- X - lambda*Y
929        and evaluate G = function(Y) (depends on the line search).
930     */
931     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
932     ynorm = 1; gnorm = fnorm;
933     ierr = (*vi->LineSearch)(snes,vi->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
934     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
935     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
936     if (snes->domainerror) {
937       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
938       PetscFunctionReturn(0);
939     }
940     if (!lssucceed) {
941       if (++snes->numFailures >= snes->maxFailures) {
942 	PetscBool ismin;
943         snes->reason = SNES_DIVERGED_LINE_SEARCH;
944         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
945         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
946         break;
947       }
948     }
949     /* Update function and solution vectors */
950     fnorm = gnorm;
951     ierr = VecCopy(G,F);CHKERRQ(ierr);
952     ierr = VecCopy(W,X);CHKERRQ(ierr);
953     /* Monitor convergence */
954     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
955     snes->iter = i+1;
956     snes->norm = fnorm;
957     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
958     SNESLogConvHistory(snes,snes->norm,lits);
959     SNESMonitor(snes,snes->iter,snes->norm);
960     /* Test for convergence, xnorm = || X || */
961     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
962     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
963     if (snes->reason) break;
964   }
965   if (i == maxits) {
966     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
967     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
968   }
969   PetscFunctionReturn(0);
970 }
971 
972 /* -------------------------------------------------------------------------- */
973 /*
974    SNESSetUp_VI - Sets up the internal data structures for the later use
975    of the SNESVI nonlinear solver.
976 
977    Input Parameter:
978 .  snes - the SNES context
979 .  x - the solution vector
980 
981    Application Interface Routine: SNESSetUp()
982 
983    Notes:
984    For basic use of the SNES solvers, the user need not explicitly call
985    SNESSetUp(), since these actions will automatically occur during
986    the call to SNESSolve().
987  */
988 #undef __FUNCT__
989 #define __FUNCT__ "SNESSetUp_VI"
990 PetscErrorCode SNESSetUp_VI(SNES snes)
991 {
992   PetscErrorCode ierr;
993   SNES_VI        *vi = (SNES_VI*) snes->data;
994   PetscInt       i_start[3],i_end[3];
995 
996   PetscFunctionBegin;
997   if (!snes->vec_sol_update) {
998     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
999     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
1000   }
1001   if (!snes->work) {
1002     snes->nwork = 3;
1003     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
1004     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
1005   }
1006 
1007   /* If the lower and upper bound on variables are not set, set it to
1008      -Inf and Inf */
1009   if (!vi->xl && !vi->xu) {
1010     vi->usersetxbounds = PETSC_FALSE;
1011     ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr);
1012     ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr);
1013     ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr);
1014     ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr);
1015   } else {
1016     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
1017     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
1018     ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr);
1019     ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr);
1020     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]))
1021       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.");
1022   }
1023   if (snes->ops->solve == SNESSolveVI_RS) {
1024     /* Set up previous active index set for the first snes solve
1025        vi->IS_inact_prev = 0,1,2,....N */
1026     PetscInt *indices;
1027     PetscInt  i,n,rstart,rend;
1028 
1029     ierr = VecGetOwnershipRange(snes->vec_sol,&rstart,&rend);CHKERRQ(ierr);
1030     ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
1031     ierr = PetscMalloc(n*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1032     for(i=0;i < n; i++) indices[i] = rstart + i;
1033     ierr = ISCreateGeneral(((PetscObject)snes)->comm,n,indices,PETSC_OWN_POINTER,&vi->IS_inact_prev);
1034   }
1035 
1036   if (snes->ops->solve != SNESSolveVI_RS) {
1037     ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr);
1038     ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr);
1039     ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr);
1040     ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr);
1041     ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr);
1042     ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr);
1043   }
1044   PetscFunctionReturn(0);
1045 }
1046 /* -------------------------------------------------------------------------- */
1047 /*
1048    SNESDestroy_VI - Destroys the private SNES_VI context that was created
1049    with SNESCreate_VI().
1050 
1051    Input Parameter:
1052 .  snes - the SNES context
1053 
1054    Application Interface Routine: SNESDestroy()
1055  */
1056 #undef __FUNCT__
1057 #define __FUNCT__ "SNESDestroy_VI"
1058 PetscErrorCode SNESDestroy_VI(SNES snes)
1059 {
1060   SNES_VI        *vi = (SNES_VI*) snes->data;
1061   PetscErrorCode ierr;
1062 
1063   PetscFunctionBegin;
1064   if (snes->vec_sol_update) {
1065     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
1066     snes->vec_sol_update = PETSC_NULL;
1067   }
1068   if (snes->nwork) {
1069     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
1070     snes->nwork = 0;
1071     snes->work  = PETSC_NULL;
1072   }
1073   if (snes->ops->solve == SNESSolveVI_RS) {
1074     ierr = ISDestroy(vi->IS_inact_prev);
1075   }
1076 
1077   if (snes->ops->solve != SNESSolveVI_RS) {
1078     /* clear vectors */
1079     ierr = VecDestroy(vi->dpsi);CHKERRQ(ierr);
1080     ierr = VecDestroy(vi->phi); CHKERRQ(ierr);
1081     ierr = VecDestroy(vi->Da); CHKERRQ(ierr);
1082     ierr = VecDestroy(vi->Db); CHKERRQ(ierr);
1083     ierr = VecDestroy(vi->z); CHKERRQ(ierr);
1084     ierr = VecDestroy(vi->t); CHKERRQ(ierr);
1085   }
1086 
1087   if (!vi->usersetxbounds) {
1088     ierr = VecDestroy(vi->xl); CHKERRQ(ierr);
1089     ierr = VecDestroy(vi->xu); CHKERRQ(ierr);
1090   }
1091 
1092   if (vi->lsmonitor) {
1093     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
1094   }
1095   ierr = PetscFree(snes->data);CHKERRQ(ierr);
1096 
1097   /* clear composed functions */
1098   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
1099   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
1100   PetscFunctionReturn(0);
1101 }
1102 
1103 /* -------------------------------------------------------------------------- */
1104 #undef __FUNCT__
1105 #define __FUNCT__ "SNESLineSearchNo_VI"
1106 
1107 /*
1108   This routine does not actually do a line search but it takes a full newton
1109   step while ensuring that the new iterates remain within the constraints.
1110 
1111 */
1112 PetscErrorCode SNESLineSearchNo_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
1113 {
1114   PetscErrorCode ierr;
1115   SNES_VI        *vi = (SNES_VI*)snes->data;
1116   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1117 
1118   PetscFunctionBegin;
1119   *flag = PETSC_TRUE;
1120   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1121   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
1122   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1123   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1124   if (vi->postcheckstep) {
1125    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1126   }
1127   if (changed_y) {
1128     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1129     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1130   }
1131   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1132   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1133   if (!snes->domainerror) {
1134     if (snes->ops->solve == SNESSolveVI_RS) {
1135        ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1136     } else {
1137       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
1138     }
1139     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1140   }
1141   if (vi->lsmonitor) {
1142     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1143   }
1144   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1145   PetscFunctionReturn(0);
1146 }
1147 
1148 /* -------------------------------------------------------------------------- */
1149 #undef __FUNCT__
1150 #define __FUNCT__ "SNESLineSearchNoNorms_VI"
1151 
1152 /*
1153   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
1154 */
1155 PetscErrorCode SNESLineSearchNoNorms_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
1156 {
1157   PetscErrorCode ierr;
1158   SNES_VI        *vi = (SNES_VI*)snes->data;
1159   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1160 
1161   PetscFunctionBegin;
1162   *flag = PETSC_TRUE;
1163   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1164   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
1165   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1166   if (vi->postcheckstep) {
1167    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1168   }
1169   if (changed_y) {
1170     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1171     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1172   }
1173 
1174   /* don't evaluate function the last time through */
1175   if (snes->iter < snes->max_its-1) {
1176     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1177   }
1178   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1179   PetscFunctionReturn(0);
1180 }
1181 
1182 /* -------------------------------------------------------------------------- */
1183 #undef __FUNCT__
1184 #define __FUNCT__ "SNESLineSearchCubic_VI"
1185 /*
1186   This routine implements a cubic line search while doing a projection on the variable bounds
1187 */
1188 PetscErrorCode SNESLineSearchCubic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
1189 {
1190   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
1191   PetscReal      minlambda,lambda,lambdatemp;
1192 #if defined(PETSC_USE_COMPLEX)
1193   PetscScalar    cinitslope;
1194 #endif
1195   PetscErrorCode ierr;
1196   PetscInt       count;
1197   SNES_VI        *vi = (SNES_VI*)snes->data;
1198   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1199   MPI_Comm       comm;
1200 
1201   PetscFunctionBegin;
1202   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1203   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1204   *flag   = PETSC_TRUE;
1205 
1206   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
1207   if (*ynorm == 0.0) {
1208     if (vi->lsmonitor) {
1209       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
1210     }
1211     *gnorm = fnorm;
1212     ierr   = VecCopy(x,w);CHKERRQ(ierr);
1213     ierr   = VecCopy(f,g);CHKERRQ(ierr);
1214     *flag  = PETSC_FALSE;
1215     goto theend1;
1216   }
1217   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
1218     if (vi->lsmonitor) {
1219       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr);
1220     }
1221     ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
1222     *ynorm = vi->maxstep;
1223   }
1224   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
1225   minlambda = vi->minlambda/rellength;
1226   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
1227 #if defined(PETSC_USE_COMPLEX)
1228   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
1229   initslope = PetscRealPart(cinitslope);
1230 #else
1231   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
1232 #endif
1233   if (initslope > 0.0)  initslope = -initslope;
1234   if (initslope == 0.0) initslope = -1.0;
1235 
1236   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1237   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1238   if (snes->nfuncs >= snes->max_funcs) {
1239     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
1240     *flag = PETSC_FALSE;
1241     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1242     goto theend1;
1243   }
1244   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1245   if (snes->ops->solve == SNESSolveVI_RS) {
1246     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1247   } else {
1248     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1249   }
1250   if (snes->domainerror) {
1251     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1252     PetscFunctionReturn(0);
1253   }
1254   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1255   ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1256   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1257     if (vi->lsmonitor) {
1258       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1259     }
1260     goto theend1;
1261   }
1262 
1263   /* Fit points with quadratic */
1264   lambda     = 1.0;
1265   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
1266   lambdaprev = lambda;
1267   gnormprev  = *gnorm;
1268   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
1269   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
1270   else                         lambda = lambdatemp;
1271 
1272   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
1273   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1274   if (snes->nfuncs >= snes->max_funcs) {
1275     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
1276     *flag = PETSC_FALSE;
1277     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1278     goto theend1;
1279   }
1280   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1281   if (snes->ops->solve == SNESSolveVI_RS) {
1282     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1283   } else {
1284     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1285   }
1286   if (snes->domainerror) {
1287     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1288     PetscFunctionReturn(0);
1289   }
1290   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1291   if (vi->lsmonitor) {
1292     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr);
1293   }
1294   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1295     if (vi->lsmonitor) {
1296       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
1297     }
1298     goto theend1;
1299   }
1300 
1301   /* Fit points with cubic */
1302   count = 1;
1303   while (PETSC_TRUE) {
1304     if (lambda <= minlambda) {
1305       if (vi->lsmonitor) {
1306 	ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
1307 	ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,minlambda,lambda,initslope);CHKERRQ(ierr);
1308       }
1309       *flag = PETSC_FALSE;
1310       break;
1311     }
1312     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
1313     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
1314     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
1315     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
1316     d  = b*b - 3*a*initslope;
1317     if (d < 0.0) d = 0.0;
1318     if (a == 0.0) {
1319       lambdatemp = -initslope/(2.0*b);
1320     } else {
1321       lambdatemp = (-b + sqrt(d))/(3.0*a);
1322     }
1323     lambdaprev = lambda;
1324     gnormprev  = *gnorm;
1325     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
1326     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
1327     else                         lambda     = lambdatemp;
1328     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
1329     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1330     if (snes->nfuncs >= snes->max_funcs) {
1331       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
1332       ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
1333       *flag = PETSC_FALSE;
1334       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1335       break;
1336     }
1337     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1338     if (snes->ops->solve == SNESSolveVI_RS) {
1339       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1340     } else {
1341       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1342     }
1343     if (snes->domainerror) {
1344       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1345       PetscFunctionReturn(0);
1346     }
1347     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1348     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */
1349       if (vi->lsmonitor) {
1350 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
1351       }
1352       break;
1353     } else {
1354       if (vi->lsmonitor) {
1355         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
1356       }
1357     }
1358     count++;
1359   }
1360   theend1:
1361   /* Optional user-defined check for line search step validity */
1362   if (vi->postcheckstep && *flag) {
1363     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1364     if (changed_y) {
1365       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1366       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1367     }
1368     if (changed_y || changed_w) { /* recompute the function if the step has changed */
1369       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1370       if (snes->ops->solve == SNESSolveVI_RS) {
1371         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1372       } else {
1373         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1374       }
1375       if (snes->domainerror) {
1376         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1377         PetscFunctionReturn(0);
1378       }
1379       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1380       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
1381       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
1382     }
1383   }
1384   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1385   PetscFunctionReturn(0);
1386 }
1387 
1388 #undef __FUNCT__
1389 #define __FUNCT__ "SNESLineSearchQuadratic_VI"
1390 /*
1391   This routine does a quadratic line search while keeping the iterates within the variable bounds
1392 */
1393 PetscErrorCode SNESLineSearchQuadratic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
1394 {
1395   /*
1396      Note that for line search purposes we work with with the related
1397      minimization problem:
1398         min  z(x):  R^n -> R,
1399      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
1400    */
1401   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
1402 #if defined(PETSC_USE_COMPLEX)
1403   PetscScalar    cinitslope;
1404 #endif
1405   PetscErrorCode ierr;
1406   PetscInt       count;
1407   SNES_VI        *vi = (SNES_VI*)snes->data;
1408   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1409 
1410   PetscFunctionBegin;
1411   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1412   *flag   = PETSC_TRUE;
1413 
1414   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
1415   if (*ynorm == 0.0) {
1416     if (vi->lsmonitor) {
1417       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
1418     }
1419     *gnorm = fnorm;
1420     ierr   = VecCopy(x,w);CHKERRQ(ierr);
1421     ierr   = VecCopy(f,g);CHKERRQ(ierr);
1422     *flag  = PETSC_FALSE;
1423     goto theend2;
1424   }
1425   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
1426     ierr   = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
1427     *ynorm = vi->maxstep;
1428   }
1429   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
1430   minlambda = vi->minlambda/rellength;
1431   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
1432 #if defined(PETSC_USE_COMPLEX)
1433   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
1434   initslope = PetscRealPart(cinitslope);
1435 #else
1436   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
1437 #endif
1438   if (initslope > 0.0)  initslope = -initslope;
1439   if (initslope == 0.0) initslope = -1.0;
1440   ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr);
1441 
1442   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1443   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1444   if (snes->nfuncs >= snes->max_funcs) {
1445     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
1446     *flag = PETSC_FALSE;
1447     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1448     goto theend2;
1449   }
1450   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1451   if (snes->ops->solve == SNESSolveVI_RS) {
1452     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1453   } else {
1454     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1455   }
1456   if (snes->domainerror) {
1457     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1458     PetscFunctionReturn(0);
1459   }
1460   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1461   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1462     if (vi->lsmonitor) {
1463       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1464     }
1465     goto theend2;
1466   }
1467 
1468   /* Fit points with quadratic */
1469   lambda = 1.0;
1470   count = 1;
1471   while (PETSC_TRUE) {
1472     if (lambda <= minlambda) { /* bad luck; use full step */
1473       if (vi->lsmonitor) {
1474         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
1475         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
1476       }
1477       ierr = VecCopy(x,w);CHKERRQ(ierr);
1478       *flag = PETSC_FALSE;
1479       break;
1480     }
1481     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
1482     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
1483     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
1484     else                         lambda     = lambdatemp;
1485 
1486     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
1487     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1488     if (snes->nfuncs >= snes->max_funcs) {
1489       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
1490       ierr  = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
1491       *flag = PETSC_FALSE;
1492       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1493       break;
1494     }
1495     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1496     if (snes->domainerror) {
1497       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1498       PetscFunctionReturn(0);
1499     }
1500     if (snes->ops->solve == SNESSolveVI_RS) {
1501       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1502     } else {
1503       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1504     }
1505     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1506     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1507       if (vi->lsmonitor) {
1508         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
1509       }
1510       break;
1511     }
1512     count++;
1513   }
1514   theend2:
1515   /* Optional user-defined check for line search step validity */
1516   if (vi->postcheckstep) {
1517     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1518     if (changed_y) {
1519       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1520       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1521     }
1522     if (changed_y || changed_w) { /* recompute the function if the step has changed */
1523       ierr = SNESComputeFunction(snes,w,g);
1524       if (snes->domainerror) {
1525         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1526         PetscFunctionReturn(0);
1527       }
1528       if (snes->ops->solve == SNESSolveVI_RS) {
1529         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1530       } else {
1531         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1532       }
1533 
1534       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
1535       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
1536       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1537     }
1538   }
1539   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1540   PetscFunctionReturn(0);
1541 }
1542 
1543 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/
1544 /* -------------------------------------------------------------------------- */
1545 EXTERN_C_BEGIN
1546 #undef __FUNCT__
1547 #define __FUNCT__ "SNESLineSearchSet_VI"
1548 PetscErrorCode  SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx)
1549 {
1550   PetscFunctionBegin;
1551   ((SNES_VI *)(snes->data))->LineSearch = func;
1552   ((SNES_VI *)(snes->data))->lsP        = lsctx;
1553   PetscFunctionReturn(0);
1554 }
1555 EXTERN_C_END
1556 
1557 /* -------------------------------------------------------------------------- */
1558 EXTERN_C_BEGIN
1559 #undef __FUNCT__
1560 #define __FUNCT__ "SNESLineSearchSetMonitor_VI"
1561 PetscErrorCode  SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg)
1562 {
1563   SNES_VI        *vi = (SNES_VI*)snes->data;
1564   PetscErrorCode ierr;
1565 
1566   PetscFunctionBegin;
1567   if (flg && !vi->lsmonitor) {
1568     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr);
1569   } else if (!flg && vi->lsmonitor) {
1570     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
1571   }
1572   PetscFunctionReturn(0);
1573 }
1574 EXTERN_C_END
1575 
1576 /*
1577    SNESView_VI - Prints info from the SNESVI data structure.
1578 
1579    Input Parameters:
1580 .  SNES - the SNES context
1581 .  viewer - visualization context
1582 
1583    Application Interface Routine: SNESView()
1584 */
1585 #undef __FUNCT__
1586 #define __FUNCT__ "SNESView_VI"
1587 static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer)
1588 {
1589   SNES_VI        *vi = (SNES_VI *)snes->data;
1590   const char     *cstr,*tstr;
1591   PetscErrorCode ierr;
1592   PetscBool     iascii;
1593 
1594   PetscFunctionBegin;
1595   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1596   if (iascii) {
1597     if (vi->LineSearch == SNESLineSearchNo_VI)             cstr = "SNESLineSearchNo";
1598     else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic";
1599     else if (vi->LineSearch == SNESLineSearchCubic_VI)     cstr = "SNESLineSearchCubic";
1600     else                                                             cstr = "unknown";
1601     if (snes->ops->solve == SNESSolveVI_SS)      tstr = "Semismooth";
1602     else if (snes->ops->solve == SNESSolveVI_RS) tstr = "Reduced Space";
1603     else                                         tstr = "unknown";
1604     ierr = PetscViewerASCIIPrintf(viewer,"  VI algorithm: %s\n",tstr);CHKERRQ(ierr);
1605     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1606     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr);
1607   } else {
1608     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name);
1609   }
1610   PetscFunctionReturn(0);
1611 }
1612 
1613 #undef __FUNCT__
1614 #define __FUNCT__ "SNESVISetVariableBounds"
1615 /*@
1616    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
1617 
1618    Input Parameters:
1619 .  snes - the SNES context.
1620 .  xl   - lower bound.
1621 .  xu   - upper bound.
1622 
1623    Notes:
1624    If this routine is not called then the lower and upper bounds are set to
1625    PETSC_VI_INF and PETSC_VI_NINF respectively during SNESSetUp().
1626 
1627 @*/
1628 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
1629 {
1630   SNES_VI  *vi = (SNES_VI*)snes->data;
1631 
1632   PetscFunctionBegin;
1633   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1634   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
1635   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
1636   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1637   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);
1638   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);
1639 
1640   vi->usersetxbounds = PETSC_TRUE;
1641   vi->xl = xl;
1642   vi->xu = xu;
1643   PetscFunctionReturn(0);
1644 }
1645 
1646 /* -------------------------------------------------------------------------- */
1647 /*
1648    SNESSetFromOptions_VI - Sets various parameters for the SNESVI method.
1649 
1650    Input Parameter:
1651 .  snes - the SNES context
1652 
1653    Application Interface Routine: SNESSetFromOptions()
1654 */
1655 #undef __FUNCT__
1656 #define __FUNCT__ "SNESSetFromOptions_VI"
1657 static PetscErrorCode SNESSetFromOptions_VI(SNES snes)
1658 {
1659   SNES_VI        *vi = (SNES_VI *)snes->data;
1660   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1661   const char     *vies[] = {"ss","rs"};
1662   PetscErrorCode ierr;
1663   PetscInt       indx;
1664   PetscBool      flg,set,flg2;
1665 
1666   PetscFunctionBegin;
1667   ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr);
1668   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
1669   if (flg) {
1670     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
1671   }
1672   ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr);
1673   ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr);
1674   ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr);
1675   ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr);
1676   ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
1677   if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
1678   ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr);
1679   if (flg2) {
1680     switch (indx) {
1681     case 0:
1682       snes->ops->solve = SNESSolveVI_SS;
1683       break;
1684     case 1:
1685       snes->ops->solve = SNESSolveVI_RS;
1686       break;
1687     }
1688   }
1689   ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"basic",&indx,&flg);CHKERRQ(ierr);
1690   if (flg) {
1691     switch (indx) {
1692     case 0:
1693       ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
1694       break;
1695     case 1:
1696       ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
1697       break;
1698     case 2:
1699       ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
1700       break;
1701     case 3:
1702       ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
1703       break;
1704     }
1705   }
1706   ierr = PetscOptionsTail();CHKERRQ(ierr);
1707   PetscFunctionReturn(0);
1708 }
1709 /* -------------------------------------------------------------------------- */
1710 /*MC
1711       SNESVI - Semismooth newton method based nonlinear solver that uses a line search
1712 
1713    Options Database:
1714 +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1715 .   -snes_ls_alpha <alpha> - Sets alpha
1716 .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
1717 .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
1718 -   -snes_ls_monitor - print information about progress of line searches
1719 
1720 
1721    Level: beginner
1722 
1723 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
1724            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
1725            SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
1726 
1727 M*/
1728 EXTERN_C_BEGIN
1729 #undef __FUNCT__
1730 #define __FUNCT__ "SNESCreate_VI"
1731 PetscErrorCode  SNESCreate_VI(SNES snes)
1732 {
1733   PetscErrorCode ierr;
1734   SNES_VI      *vi;
1735 
1736   PetscFunctionBegin;
1737   snes->ops->setup	     = SNESSetUp_VI;
1738   snes->ops->solve	     = SNESSolveVI_SS;
1739   snes->ops->destroy	     = SNESDestroy_VI;
1740   snes->ops->setfromoptions  = SNESSetFromOptions_VI;
1741   snes->ops->view            = SNESView_VI;
1742   snes->ops->converged       = SNESDefaultConverged_VI;
1743 
1744   ierr                   = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr);
1745   snes->data    	 = (void*)vi;
1746   vi->alpha		 = 1.e-4;
1747   vi->maxstep		 = 1.e8;
1748   vi->minlambda         = 1.e-12;
1749   vi->LineSearch        = SNESLineSearchNo_VI;
1750   vi->lsP               = PETSC_NULL;
1751   vi->postcheckstep     = PETSC_NULL;
1752   vi->postcheck         = PETSC_NULL;
1753   vi->precheckstep      = PETSC_NULL;
1754   vi->precheck          = PETSC_NULL;
1755   vi->const_tol         =  2.2204460492503131e-16;
1756 
1757   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr);
1758   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr);
1759 
1760   PetscFunctionReturn(0);
1761 }
1762 EXTERN_C_END
1763