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