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