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