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