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