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