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