xref: /petsc/src/snes/impls/vi/vi.c (revision 6ce1633cb736e3bd2a11b0bc146401a5bd4cb96c)
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,Nis_inact_prev;
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   KSP                snesksp;
763   Mat                Amat;
764 
765   PetscFunctionBegin;
766   snes->numFailures            = 0;
767   snes->numLinearSolveFailures = 0;
768   snes->reason                 = SNES_CONVERGED_ITERATING;
769 
770   maxits	= snes->max_its;	/* maximum number of iterations */
771   X		= snes->vec_sol;	/* solution vector */
772   F		= snes->vec_func;	/* residual vector */
773   Y		= snes->work[0];	/* work vectors */
774   G		= snes->work[1];
775   W		= snes->work[2];
776 
777   /* Get the inactive set size from the previous snes solve */
778   ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr);
779   ierr = KSPGetOperators(snesksp,&Amat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
780   ierr = MatGetSize(Amat,&Nis_inact_prev,PETSC_NULL);CHKERRQ(ierr);
781 
782   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
783   snes->iter = 0;
784   snes->norm = 0.0;
785   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
786 
787   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
788   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
789   if (snes->domainerror) {
790     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
791     PetscFunctionReturn(0);
792   }
793   ierr = SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);CHKERRQ(ierr);
794   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
795   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
796   if PetscIsInfOrNanReal(fnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
797 
798   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
799   snes->norm = fnorm;
800   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
801   SNESLogConvHistory(snes,fnorm,0);
802   SNESMonitor(snes,0,fnorm);
803 
804   /* set parameter for default relative tolerance convergence test */
805   snes->ttol = fnorm*snes->rtol;
806   /* test convergence */
807   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
808   if (snes->reason) PetscFunctionReturn(0);
809 
810   for (i=0; i<maxits; i++) {
811 
812     IS         IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
813     VecScatter scat_act,scat_inact;
814     PetscInt   nis_act,nis_inact,Nis_inact;
815     Vec        Y_act,Y_inact,F_inact;
816     Mat        jac_inact_inact,prejac_inact_inact;
817     IS         keptrows;
818 
819     /* Call general purpose update function */
820     if (snes->ops->update) {
821       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
822     }
823     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
824 
825     /* Create active and inactive index sets */
826     ierr = SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&IS_inact);CHKERRQ(ierr);
827 
828     /* Create inactive set submatrices */
829     ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
830     ierr = MatFindNonzeroRows(jac_inact_inact,&keptrows);CHKERRQ(ierr);
831     if (keptrows) {
832       PetscInt       cnt,*nrows,k;
833       const PetscInt *krows,*inact;
834 
835       ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr);
836       ierr = ISDestroy(IS_act);CHKERRQ(ierr);
837 
838       ierr = ISGetLocalSize(keptrows,&cnt);CHKERRQ(ierr);
839       ierr = ISGetIndices(keptrows,&krows);CHKERRQ(ierr);
840       ierr = ISGetIndices(IS_inact,&inact);CHKERRQ(ierr);
841       ierr = PetscMalloc(cnt*sizeof(PetscInt),&nrows);CHKERRQ(ierr);
842       for (k=0; k<cnt; k++) {
843         nrows[k] = inact[krows[k]];
844       }
845       ierr = ISRestoreIndices(keptrows,&krows);CHKERRQ(ierr);
846       ierr = ISRestoreIndices(IS_inact,&inact);CHKERRQ(ierr);
847       ierr = ISDestroy(keptrows);CHKERRQ(ierr);
848       ierr = ISDestroy(IS_inact);CHKERRQ(ierr);
849 
850       ierr = ISCreateGeneral(PETSC_COMM_SELF,cnt,nrows,PETSC_OWN_POINTER,&IS_inact);CHKERRQ(ierr);
851       ierr = ISComplement(IS_inact,0,F->map->n,&IS_act);CHKERRQ(ierr);
852       ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
853     }
854 
855     /* Get sizes of active and inactive sets */
856     ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr);
857     ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr);
858     ierr = ISGetSize(IS_inact,&Nis_inact);CHKERRQ(ierr);
859 
860     /* Create active and inactive set vectors */
861     ierr = SNESVICreateSubVectors(snes,nis_inact,&F_inact);CHKERRQ(ierr);
862     ierr = SNESVICreateSubVectors(snes,nis_act,&Y_act);CHKERRQ(ierr);
863     ierr = SNESVICreateSubVectors(snes,nis_inact,&Y_inact);CHKERRQ(ierr);
864 
865     /* Create scatter contexts */
866     ierr = VecScatterCreate(Y,IS_act,Y_act,PETSC_NULL,&scat_act);CHKERRQ(ierr);
867     ierr = VecScatterCreate(Y,IS_inact,Y_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr);
868 
869     /* Do a vec scatter to active and inactive set vectors */
870     ierr = VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
871     ierr = VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
872 
873     ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
874     ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
875 
876     ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
877     ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
878 
879     /* Active set direction = 0 */
880     ierr = VecSet(Y_act,0);CHKERRQ(ierr);
881     if (snes->jacobian != snes->jacobian_pre) {
882       ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr);
883     } else prejac_inact_inact = jac_inact_inact;
884 
885     if (Nis_inact != Nis_inact_prev) {
886       ierr = SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);CHKERRQ(ierr);
887     }
888 
889     ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr);
890     ierr = SNES_KSPSolve(snes,snes->ksp,F_inact,Y_inact);CHKERRQ(ierr);
891     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
892     if (kspreason < 0) {
893       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
894         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
895         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
896         break;
897       }
898      }
899 
900     ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
901     ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
902     ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
903     ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
904 
905     ierr = VecDestroy(F_inact);CHKERRQ(ierr);
906     ierr = VecDestroy(Y_act);CHKERRQ(ierr);
907     ierr = VecDestroy(Y_inact);CHKERRQ(ierr);
908     ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr);
909     ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr);
910     ierr = ISDestroy(IS_act);CHKERRQ(ierr);
911     ierr = ISDestroy(IS_inact);CHKERRQ(ierr);
912     ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr);
913     if (snes->jacobian != snes->jacobian_pre) {
914       ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr);
915     }
916     Nis_inact_prev = Nis_inact;
917     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
918     snes->linear_its += lits;
919     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
920     /*
921     if (vi->precheckstep) {
922       PetscBool changed_y = PETSC_FALSE;
923       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
924     }
925 
926     if (PetscLogPrintInfo){
927       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
928     }
929     */
930     /* Compute a (scaled) negative update in the line search routine:
931          Y <- X - lambda*Y
932        and evaluate G = function(Y) (depends on the line search).
933     */
934     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
935     ynorm = 1; gnorm = fnorm;
936     ierr = (*vi->LineSearch)(snes,vi->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
937     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
938     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
939     if (snes->domainerror) {
940       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
941       PetscFunctionReturn(0);
942     }
943     if (!lssucceed) {
944       if (++snes->numFailures >= snes->maxFailures) {
945 	PetscBool ismin;
946         snes->reason = SNES_DIVERGED_LINE_SEARCH;
947         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
948         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
949         break;
950       }
951     }
952     /* Update function and solution vectors */
953     fnorm = gnorm;
954     ierr = VecCopy(G,F);CHKERRQ(ierr);
955     ierr = VecCopy(W,X);CHKERRQ(ierr);
956     /* Monitor convergence */
957     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
958     snes->iter = i+1;
959     snes->norm = fnorm;
960     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
961     SNESLogConvHistory(snes,snes->norm,lits);
962     SNESMonitor(snes,snes->iter,snes->norm);
963     /* Test for convergence, xnorm = || X || */
964     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
965     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
966     if (snes->reason) break;
967   }
968   if (i == maxits) {
969     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
970     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
971   }
972   PetscFunctionReturn(0);
973 }
974 
975 /* -------------------------------------------------------------------------- */
976 /*
977    SNESSetUp_VI - Sets up the internal data structures for the later use
978    of the SNESVI nonlinear solver.
979 
980    Input Parameter:
981 .  snes - the SNES context
982 .  x - the solution vector
983 
984    Application Interface Routine: SNESSetUp()
985 
986    Notes:
987    For basic use of the SNES solvers, the user need not explicitly call
988    SNESSetUp(), since these actions will automatically occur during
989    the call to SNESSolve().
990  */
991 #undef __FUNCT__
992 #define __FUNCT__ "SNESSetUp_VI"
993 PetscErrorCode SNESSetUp_VI(SNES snes)
994 {
995   PetscErrorCode ierr;
996   SNES_VI        *vi = (SNES_VI*) snes->data;
997   PetscInt       i_start[3],i_end[3];
998 
999   PetscFunctionBegin;
1000   if (!snes->vec_sol_update) {
1001     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
1002     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
1003   }
1004   if (!snes->work) {
1005     snes->nwork = 3;
1006     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
1007     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
1008   }
1009 
1010   /* If the lower and upper bound on variables are not set, set it to
1011      -Inf and Inf */
1012   if (!vi->xl && !vi->xu) {
1013     vi->usersetxbounds = PETSC_FALSE;
1014     ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr);
1015     ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr);
1016     ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr);
1017     ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr);
1018   } else {
1019     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
1020     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
1021     ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr);
1022     ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr);
1023     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]))
1024       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.");
1025   }
1026 
1027   if (snes->ops->solve != SNESSolveVI_RS) {
1028     ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr);
1029     ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr);
1030     ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr);
1031     ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr);
1032     ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr);
1033     ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr);
1034 
1035   }
1036   PetscFunctionReturn(0);
1037 }
1038 /* -------------------------------------------------------------------------- */
1039 /*
1040    SNESDestroy_VI - Destroys the private SNES_VI context that was created
1041    with SNESCreate_VI().
1042 
1043    Input Parameter:
1044 .  snes - the SNES context
1045 
1046    Application Interface Routine: SNESDestroy()
1047  */
1048 #undef __FUNCT__
1049 #define __FUNCT__ "SNESDestroy_VI"
1050 PetscErrorCode SNESDestroy_VI(SNES snes)
1051 {
1052   SNES_VI        *vi = (SNES_VI*) snes->data;
1053   PetscErrorCode ierr;
1054 
1055   PetscFunctionBegin;
1056   if (snes->vec_sol_update) {
1057     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
1058     snes->vec_sol_update = PETSC_NULL;
1059   }
1060   if (snes->nwork) {
1061     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
1062     snes->nwork = 0;
1063     snes->work  = PETSC_NULL;
1064   }
1065 
1066   if (snes->ops->solve != SNESSolveVI_RS) {
1067     /* clear vectors */
1068     ierr = VecDestroy(vi->dpsi);CHKERRQ(ierr);
1069     ierr = VecDestroy(vi->phi); CHKERRQ(ierr);
1070     ierr = VecDestroy(vi->Da); CHKERRQ(ierr);
1071     ierr = VecDestroy(vi->Db); CHKERRQ(ierr);
1072     ierr = VecDestroy(vi->z); CHKERRQ(ierr);
1073     ierr = VecDestroy(vi->t); CHKERRQ(ierr);
1074   }
1075 
1076   if (!vi->usersetxbounds) {
1077     ierr = VecDestroy(vi->xl); CHKERRQ(ierr);
1078     ierr = VecDestroy(vi->xu); CHKERRQ(ierr);
1079   }
1080 
1081   if (vi->lsmonitor) {
1082     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
1083   }
1084   ierr = PetscFree(snes->data);CHKERRQ(ierr);
1085 
1086   /* clear composed functions */
1087   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
1088   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
1089   PetscFunctionReturn(0);
1090 }
1091 
1092 /* -------------------------------------------------------------------------- */
1093 #undef __FUNCT__
1094 #define __FUNCT__ "SNESLineSearchNo_VI"
1095 
1096 /*
1097   This routine does not actually do a line search but it takes a full newton
1098   step while ensuring that the new iterates remain within the constraints.
1099 
1100 */
1101 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)
1102 {
1103   PetscErrorCode ierr;
1104   SNES_VI        *vi = (SNES_VI*)snes->data;
1105   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1106 
1107   PetscFunctionBegin;
1108   *flag = PETSC_TRUE;
1109   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1110   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
1111   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1112   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1113   if (vi->postcheckstep) {
1114    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1115   }
1116   if (changed_y) {
1117     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1118     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1119   }
1120   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1121   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1122   if (!snes->domainerror) {
1123     if (snes->ops->solve == SNESSolveVI_RS) {
1124        ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1125     } else {
1126       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
1127     }
1128     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1129   }
1130   if (vi->lsmonitor) {
1131     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1132   }
1133   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1134   PetscFunctionReturn(0);
1135 }
1136 
1137 /* -------------------------------------------------------------------------- */
1138 #undef __FUNCT__
1139 #define __FUNCT__ "SNESLineSearchNoNorms_VI"
1140 
1141 /*
1142   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
1143 */
1144 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)
1145 {
1146   PetscErrorCode ierr;
1147   SNES_VI        *vi = (SNES_VI*)snes->data;
1148   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1149 
1150   PetscFunctionBegin;
1151   *flag = PETSC_TRUE;
1152   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1153   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
1154   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1155   if (vi->postcheckstep) {
1156    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1157   }
1158   if (changed_y) {
1159     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1160     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1161   }
1162 
1163   /* don't evaluate function the last time through */
1164   if (snes->iter < snes->max_its-1) {
1165     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1166   }
1167   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1168   PetscFunctionReturn(0);
1169 }
1170 
1171 /* -------------------------------------------------------------------------- */
1172 #undef __FUNCT__
1173 #define __FUNCT__ "SNESLineSearchCubic_VI"
1174 /*
1175   This routine implements a cubic line search while doing a projection on the variable bounds
1176 */
1177 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)
1178 {
1179   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
1180   PetscReal      minlambda,lambda,lambdatemp;
1181 #if defined(PETSC_USE_COMPLEX)
1182   PetscScalar    cinitslope;
1183 #endif
1184   PetscErrorCode ierr;
1185   PetscInt       count;
1186   SNES_VI        *vi = (SNES_VI*)snes->data;
1187   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1188   MPI_Comm       comm;
1189 
1190   PetscFunctionBegin;
1191   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1192   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1193   *flag   = PETSC_TRUE;
1194 
1195   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
1196   if (*ynorm == 0.0) {
1197     if (vi->lsmonitor) {
1198       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
1199     }
1200     *gnorm = fnorm;
1201     ierr   = VecCopy(x,w);CHKERRQ(ierr);
1202     ierr   = VecCopy(f,g);CHKERRQ(ierr);
1203     *flag  = PETSC_FALSE;
1204     goto theend1;
1205   }
1206   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
1207     if (vi->lsmonitor) {
1208       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr);
1209     }
1210     ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
1211     *ynorm = vi->maxstep;
1212   }
1213   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
1214   minlambda = vi->minlambda/rellength;
1215   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
1216 #if defined(PETSC_USE_COMPLEX)
1217   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
1218   initslope = PetscRealPart(cinitslope);
1219 #else
1220   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
1221 #endif
1222   if (initslope > 0.0)  initslope = -initslope;
1223   if (initslope == 0.0) initslope = -1.0;
1224 
1225   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1226   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1227   if (snes->nfuncs >= snes->max_funcs) {
1228     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
1229     *flag = PETSC_FALSE;
1230     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1231     goto theend1;
1232   }
1233   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1234   if (snes->ops->solve == SNESSolveVI_RS) {
1235     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1236   } else {
1237     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1238   }
1239   if (snes->domainerror) {
1240     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1241     PetscFunctionReturn(0);
1242   }
1243   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1244   ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1245   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1246     if (vi->lsmonitor) {
1247       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1248     }
1249     goto theend1;
1250   }
1251 
1252   /* Fit points with quadratic */
1253   lambda     = 1.0;
1254   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
1255   lambdaprev = lambda;
1256   gnormprev  = *gnorm;
1257   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
1258   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
1259   else                         lambda = lambdatemp;
1260 
1261   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
1262   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1263   if (snes->nfuncs >= snes->max_funcs) {
1264     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
1265     *flag = PETSC_FALSE;
1266     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1267     goto theend1;
1268   }
1269   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1270   if (snes->ops->solve == SNESSolveVI_RS) {
1271     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1272   } else {
1273     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1274   }
1275   if (snes->domainerror) {
1276     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1277     PetscFunctionReturn(0);
1278   }
1279   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1280   if (vi->lsmonitor) {
1281     ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr);
1282   }
1283   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1284     if (vi->lsmonitor) {
1285       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
1286     }
1287     goto theend1;
1288   }
1289 
1290   /* Fit points with cubic */
1291   count = 1;
1292   while (PETSC_TRUE) {
1293     if (lambda <= minlambda) {
1294       if (vi->lsmonitor) {
1295 	ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
1296 	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);
1297       }
1298       *flag = PETSC_FALSE;
1299       break;
1300     }
1301     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
1302     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
1303     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
1304     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
1305     d  = b*b - 3*a*initslope;
1306     if (d < 0.0) d = 0.0;
1307     if (a == 0.0) {
1308       lambdatemp = -initslope/(2.0*b);
1309     } else {
1310       lambdatemp = (-b + sqrt(d))/(3.0*a);
1311     }
1312     lambdaprev = lambda;
1313     gnormprev  = *gnorm;
1314     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
1315     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
1316     else                         lambda     = lambdatemp;
1317     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
1318     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1319     if (snes->nfuncs >= snes->max_funcs) {
1320       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
1321       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);
1322       *flag = PETSC_FALSE;
1323       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1324       break;
1325     }
1326     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1327     if (snes->ops->solve == SNESSolveVI_RS) {
1328       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1329     } else {
1330       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1331     }
1332     if (snes->domainerror) {
1333       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1334       PetscFunctionReturn(0);
1335     }
1336     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1337     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */
1338       if (vi->lsmonitor) {
1339 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
1340       }
1341       break;
1342     } else {
1343       if (vi->lsmonitor) {
1344         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
1345       }
1346     }
1347     count++;
1348   }
1349   theend1:
1350   /* Optional user-defined check for line search step validity */
1351   if (vi->postcheckstep && *flag) {
1352     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1353     if (changed_y) {
1354       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1355       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1356     }
1357     if (changed_y || changed_w) { /* recompute the function if the step has changed */
1358       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1359       if (snes->ops->solve == SNESSolveVI_RS) {
1360         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1361       } else {
1362         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1363       }
1364       if (snes->domainerror) {
1365         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1366         PetscFunctionReturn(0);
1367       }
1368       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1369       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
1370       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
1371     }
1372   }
1373   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1374   PetscFunctionReturn(0);
1375 }
1376 
1377 #undef __FUNCT__
1378 #define __FUNCT__ "SNESLineSearchQuadratic_VI"
1379 /*
1380   This routine does a quadratic line search while keeping the iterates within the variable bounds
1381 */
1382 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)
1383 {
1384   /*
1385      Note that for line search purposes we work with with the related
1386      minimization problem:
1387         min  z(x):  R^n -> R,
1388      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
1389    */
1390   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
1391 #if defined(PETSC_USE_COMPLEX)
1392   PetscScalar    cinitslope;
1393 #endif
1394   PetscErrorCode ierr;
1395   PetscInt       count;
1396   SNES_VI        *vi = (SNES_VI*)snes->data;
1397   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1398 
1399   PetscFunctionBegin;
1400   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1401   *flag   = PETSC_TRUE;
1402 
1403   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
1404   if (*ynorm == 0.0) {
1405     if (vi->lsmonitor) {
1406       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
1407     }
1408     *gnorm = fnorm;
1409     ierr   = VecCopy(x,w);CHKERRQ(ierr);
1410     ierr   = VecCopy(f,g);CHKERRQ(ierr);
1411     *flag  = PETSC_FALSE;
1412     goto theend2;
1413   }
1414   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
1415     ierr   = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
1416     *ynorm = vi->maxstep;
1417   }
1418   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
1419   minlambda = vi->minlambda/rellength;
1420   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
1421 #if defined(PETSC_USE_COMPLEX)
1422   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
1423   initslope = PetscRealPart(cinitslope);
1424 #else
1425   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
1426 #endif
1427   if (initslope > 0.0)  initslope = -initslope;
1428   if (initslope == 0.0) initslope = -1.0;
1429   ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr);
1430 
1431   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1432   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1433   if (snes->nfuncs >= snes->max_funcs) {
1434     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
1435     *flag = PETSC_FALSE;
1436     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1437     goto theend2;
1438   }
1439   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1440   if (snes->ops->solve == SNESSolveVI_RS) {
1441     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1442   } else {
1443     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1444   }
1445   if (snes->domainerror) {
1446     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1447     PetscFunctionReturn(0);
1448   }
1449   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1450   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1451     if (vi->lsmonitor) {
1452       ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1453     }
1454     goto theend2;
1455   }
1456 
1457   /* Fit points with quadratic */
1458   lambda = 1.0;
1459   count = 1;
1460   while (PETSC_TRUE) {
1461     if (lambda <= minlambda) { /* bad luck; use full step */
1462       if (vi->lsmonitor) {
1463         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
1464         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);
1465       }
1466       ierr = VecCopy(x,w);CHKERRQ(ierr);
1467       *flag = PETSC_FALSE;
1468       break;
1469     }
1470     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
1471     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
1472     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
1473     else                         lambda     = lambdatemp;
1474 
1475     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
1476     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1477     if (snes->nfuncs >= snes->max_funcs) {
1478       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
1479       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);
1480       *flag = PETSC_FALSE;
1481       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1482       break;
1483     }
1484     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1485     if (snes->domainerror) {
1486       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1487       PetscFunctionReturn(0);
1488     }
1489     if (snes->ops->solve == SNESSolveVI_RS) {
1490       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1491     } else {
1492       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1493     }
1494     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1495     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
1496       if (vi->lsmonitor) {
1497         ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"    Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
1498       }
1499       break;
1500     }
1501     count++;
1502   }
1503   theend2:
1504   /* Optional user-defined check for line search step validity */
1505   if (vi->postcheckstep) {
1506     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1507     if (changed_y) {
1508       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1509       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1510     }
1511     if (changed_y || changed_w) { /* recompute the function if the step has changed */
1512       ierr = SNESComputeFunction(snes,w,g);
1513       if (snes->domainerror) {
1514         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1515         PetscFunctionReturn(0);
1516       }
1517       if (snes->ops->solve == SNESSolveVI_RS) {
1518         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1519       } else {
1520         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1521       }
1522 
1523       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
1524       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
1525       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1526     }
1527   }
1528   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1529   PetscFunctionReturn(0);
1530 }
1531 
1532 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*/
1533 /* -------------------------------------------------------------------------- */
1534 EXTERN_C_BEGIN
1535 #undef __FUNCT__
1536 #define __FUNCT__ "SNESLineSearchSet_VI"
1537 PetscErrorCode  SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx)
1538 {
1539   PetscFunctionBegin;
1540   ((SNES_VI *)(snes->data))->LineSearch = func;
1541   ((SNES_VI *)(snes->data))->lsP        = lsctx;
1542   PetscFunctionReturn(0);
1543 }
1544 EXTERN_C_END
1545 
1546 /* -------------------------------------------------------------------------- */
1547 EXTERN_C_BEGIN
1548 #undef __FUNCT__
1549 #define __FUNCT__ "SNESLineSearchSetMonitor_VI"
1550 PetscErrorCode  SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg)
1551 {
1552   SNES_VI        *vi = (SNES_VI*)snes->data;
1553   PetscErrorCode ierr;
1554 
1555   PetscFunctionBegin;
1556   if (flg && !vi->lsmonitor) {
1557     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr);
1558   } else if (!flg && vi->lsmonitor) {
1559     ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr);
1560   }
1561   PetscFunctionReturn(0);
1562 }
1563 EXTERN_C_END
1564 
1565 /*
1566    SNESView_VI - Prints info from the SNESVI data structure.
1567 
1568    Input Parameters:
1569 .  SNES - the SNES context
1570 .  viewer - visualization context
1571 
1572    Application Interface Routine: SNESView()
1573 */
1574 #undef __FUNCT__
1575 #define __FUNCT__ "SNESView_VI"
1576 static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer)
1577 {
1578   SNES_VI        *vi = (SNES_VI *)snes->data;
1579   const char     *cstr,*tstr;
1580   PetscErrorCode ierr;
1581   PetscBool     iascii;
1582 
1583   PetscFunctionBegin;
1584   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1585   if (iascii) {
1586     if (vi->LineSearch == SNESLineSearchNo_VI)             cstr = "SNESLineSearchNo";
1587     else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic";
1588     else if (vi->LineSearch == SNESLineSearchCubic_VI)     cstr = "SNESLineSearchCubic";
1589     else                                                             cstr = "unknown";
1590     if (snes->ops->solve == SNESSolveVI_SS)      tstr = "Semismooth";
1591     else if (snes->ops->solve == SNESSolveVI_RS) tstr = "Reduced Space";
1592     else                                         tstr = "unknown";
1593     ierr = PetscViewerASCIIPrintf(viewer,"  VI algorithm: %s\n",tstr);CHKERRQ(ierr);
1594     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1595     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr);
1596   } else {
1597     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name);
1598   }
1599   PetscFunctionReturn(0);
1600 }
1601 
1602 #undef __FUNCT__
1603 #define __FUNCT__ "SNESVISetVariableBounds"
1604 /*@
1605    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
1606 
1607    Input Parameters:
1608 .  snes - the SNES context.
1609 .  xl   - lower bound.
1610 .  xu   - upper bound.
1611 
1612    Notes:
1613    If this routine is not called then the lower and upper bounds are set to
1614    PETSC_VI_INF and PETSC_VI_NINF respectively during SNESSetUp().
1615 
1616 @*/
1617 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
1618 {
1619   SNES_VI  *vi = (SNES_VI*)snes->data;
1620 
1621   PetscFunctionBegin;
1622   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1623   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
1624   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
1625   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1626   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);
1627   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);
1628 
1629   vi->usersetxbounds = PETSC_TRUE;
1630   vi->xl = xl;
1631   vi->xu = xu;
1632   PetscFunctionReturn(0);
1633 }
1634 
1635 /* -------------------------------------------------------------------------- */
1636 /*
1637    SNESSetFromOptions_VI - Sets various parameters for the SNESVI method.
1638 
1639    Input Parameter:
1640 .  snes - the SNES context
1641 
1642    Application Interface Routine: SNESSetFromOptions()
1643 */
1644 #undef __FUNCT__
1645 #define __FUNCT__ "SNESSetFromOptions_VI"
1646 static PetscErrorCode SNESSetFromOptions_VI(SNES snes)
1647 {
1648   SNES_VI        *vi = (SNES_VI *)snes->data;
1649   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1650   const char     *vies[] = {"ss","rs"};
1651   PetscErrorCode ierr;
1652   PetscInt       indx;
1653   PetscBool      flg,set,flg2;
1654 
1655   PetscFunctionBegin;
1656   ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr);
1657   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
1658   if (flg) {
1659     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
1660   }
1661   ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr);
1662   ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr);
1663   ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr);
1664   ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr);
1665   ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
1666   if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
1667   ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr);
1668   if (flg2) {
1669     switch (indx) {
1670     case 0:
1671       snes->ops->solve = SNESSolveVI_SS;
1672       break;
1673     case 1:
1674       snes->ops->solve = SNESSolveVI_RS;
1675       break;
1676     }
1677   }
1678   ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"basic",&indx,&flg);CHKERRQ(ierr);
1679   if (flg) {
1680     switch (indx) {
1681     case 0:
1682       ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
1683       break;
1684     case 1:
1685       ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
1686       break;
1687     case 2:
1688       ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
1689       break;
1690     case 3:
1691       ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
1692       break;
1693     }
1694   }
1695   ierr = PetscOptionsTail();CHKERRQ(ierr);
1696   PetscFunctionReturn(0);
1697 }
1698 /* -------------------------------------------------------------------------- */
1699 /*MC
1700       SNESVI - Semismooth newton method based nonlinear solver that uses a line search
1701 
1702    Options Database:
1703 +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1704 .   -snes_ls_alpha <alpha> - Sets alpha
1705 .   -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)
1706 .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
1707 -   -snes_ls_monitor - print information about progress of line searches
1708 
1709 
1710    Level: beginner
1711 
1712 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
1713            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
1714            SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
1715 
1716 M*/
1717 EXTERN_C_BEGIN
1718 #undef __FUNCT__
1719 #define __FUNCT__ "SNESCreate_VI"
1720 PetscErrorCode  SNESCreate_VI(SNES snes)
1721 {
1722   PetscErrorCode ierr;
1723   SNES_VI      *vi;
1724 
1725   PetscFunctionBegin;
1726   snes->ops->setup	     = SNESSetUp_VI;
1727   snes->ops->solve	     = SNESSolveVI_SS;
1728   snes->ops->destroy	     = SNESDestroy_VI;
1729   snes->ops->setfromoptions  = SNESSetFromOptions_VI;
1730   snes->ops->view            = SNESView_VI;
1731   snes->ops->converged       = SNESDefaultConverged_VI;
1732 
1733   ierr                   = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr);
1734   snes->data    	 = (void*)vi;
1735   vi->alpha		 = 1.e-4;
1736   vi->maxstep		 = 1.e8;
1737   vi->minlambda         = 1.e-12;
1738   vi->LineSearch        = SNESLineSearchNo_VI;
1739   vi->lsP               = PETSC_NULL;
1740   vi->postcheckstep     = PETSC_NULL;
1741   vi->postcheck         = PETSC_NULL;
1742   vi->precheckstep      = PETSC_NULL;
1743   vi->precheck          = PETSC_NULL;
1744   vi->const_tol         =  2.2204460492503131e-16;
1745 
1746   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr);
1747   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr);
1748 
1749   PetscFunctionReturn(0);
1750 }
1751 EXTERN_C_END
1752