xref: /petsc/src/snes/impls/vi/vi.c (revision 3c0bb4a45facb5092c160ad2b57affeed0b33558)
1 
2 #include <../src/snes/impls/vi/viimpl.h> /*I "petscsnes.h" I*/
3 #include <../include/private/kspimpl.h>
4 #include <../include/private/matimpl.h>
5 #include <../include/private/dmimpl.h>
6 
7 #undef __FUNCT__
8 #define __FUNCT__ "SNESVISetComputeVariableBounds"
9 /*@C
10    SNESVISetComputeVariableBounds - Sets a function  that is called to compute the variable bounds
11 
12    Input parameter
13 +  snes - the SNES context
14 -  compute - computes the bounds
15 
16 
17 @*/
18 PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
19 {
20   PetscErrorCode   ierr;
21   SNES_VI          *vi;
22 
23   PetscFunctionBegin;
24   ierr = SNESSetType(snes,SNESVI);CHKERRQ(ierr);
25   vi = (SNES_VI*)snes->data;
26   vi->computevariablebounds = compute;
27   PetscFunctionReturn(0);
28 }
29 
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "SNESVIGetInActiveSetIS"
33 /*
34    SNESVIGetInActiveSetIS - Gets the global indices for the bogus inactive set variables
35 
36    Input parameter
37 .  snes - the SNES context
38 .  X    - the snes solution vector
39 
40    Output parameter
41 .  ISact - active set index set
42 
43  */
44 PetscErrorCode SNESVIGetInActiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS* inact)
45 {
46   PetscErrorCode   ierr;
47   const PetscScalar *x,*xl,*xu,*f;
48   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
49 
50   PetscFunctionBegin;
51   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
52   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
53   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
54   ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr);
55   ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr);
56   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
57   /* Compute inactive set size */
58   for (i=0; i < nlocal;i++) {
59     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) nloc_isact++;
60   }
61 
62   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
63 
64   /* Set inactive set indices */
65   for(i=0; i < nlocal; i++) {
66     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) idx_act[i1++] = ilow+i;
67   }
68 
69    /* Create inactive set IS */
70   ierr = ISCreateGeneral(((PetscObject)upper)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr);
71 
72   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
73   ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr);
74   ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr);
75   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 /*
80     Provides a wrapper to a DM to allow it to be used to generated the interpolation/restriction from the DM for the smaller matrices and vectors
81   defined by the reduced space method.
82 
83     Simple calls the regular DM interpolation and restricts it to operation on the variables not associated with active constraints.
84 
85 */
86 typedef struct {
87   PetscInt       n;                                        /* size of vectors in the reduced DM space */
88   IS             inactive;
89   Vec            upper,lower,values,F;                    /* upper and lower bounds of all variables on this level, the values and the function values */
90   PetscErrorCode (*getinterpolation)(DM,DM,Mat*,Vec*);    /* DM's original routines */
91   PetscErrorCode (*coarsen)(DM, MPI_Comm, DM*);
92   PetscErrorCode (*createglobalvector)(DM,Vec*);
93   DM             dm;                                      /* when destroying this object we need to reset the above function into the base DM */
94 } DM_SNESVI;
95 
96 #undef __FUNCT__
97 #define __FUNCT__ "DMCreateGlobalVector_SNESVI"
98 /*
99      DMCreateGlobalVector_SNESVI - Creates global vector of the size of the reduced space
100 
101 */
102 PetscErrorCode  DMCreateGlobalVector_SNESVI(DM dm,Vec *vec)
103 {
104   PetscErrorCode          ierr;
105   PetscContainer          isnes;
106   DM_SNESVI               *dmsnesvi;
107 
108   PetscFunctionBegin;
109   ierr = PetscObjectQuery((PetscObject)dm,"VI",(PetscObject *)&isnes);CHKERRQ(ierr);
110   if (!isnes) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_PLIB,"Composed SNES is missing");
111   ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi);CHKERRQ(ierr);
112   ierr = VecCreateMPI(((PetscObject)dm)->comm,dmsnesvi->n,PETSC_DETERMINE,vec);CHKERRQ(ierr);
113   PetscFunctionReturn(0);
114 }
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "DMGetInterpolation_SNESVI"
118 /*
119      DMGetInterpolation_SNESVI - Modifieds the interpolation obtained from the DM by removing all rows and columns associated with active constraints.
120 
121 */
122 PetscErrorCode  DMGetInterpolation_SNESVI(DM dm1,DM dm2,Mat *mat,Vec *vec)
123 {
124   PetscErrorCode          ierr;
125   PetscContainer          isnes;
126   DM_SNESVI               *dmsnesvi1,*dmsnesvi2;
127   Mat                     interp;
128 
129   PetscFunctionBegin;
130   ierr = PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject *)&isnes);CHKERRQ(ierr);
131   if (!isnes) SETERRQ(((PetscObject)dm1)->comm,PETSC_ERR_PLIB,"Composed VI data structure is missing");
132   ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);CHKERRQ(ierr);
133   ierr = PetscObjectQuery((PetscObject)dm2,"VI",(PetscObject *)&isnes);CHKERRQ(ierr);
134   if (!isnes) SETERRQ(((PetscObject)dm2)->comm,PETSC_ERR_PLIB,"Composed VI data structure is missing");
135   ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi2);CHKERRQ(ierr);
136 
137   ierr = (*dmsnesvi1->getinterpolation)(dm1,dm2,&interp,PETSC_NULL);CHKERRQ(ierr);
138   ierr = MatGetSubMatrix(interp,dmsnesvi2->inactive,dmsnesvi1->inactive,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr);
139   ierr = MatDestroy(&interp);CHKERRQ(ierr);
140   *vec = 0;
141   PetscFunctionReturn(0);
142 }
143 
144 extern PetscErrorCode  DMSetVI(DM,Vec,Vec,Vec,Vec,IS);
145 
146 #undef __FUNCT__
147 #define __FUNCT__ "DMCoarsen_SNESVI"
148 /*
149      DMCoarsen_SNESVI - Computes the regular coarsened DM then computes additional information about its inactive set
150 
151 */
152 PetscErrorCode  DMCoarsen_SNESVI(DM dm1,MPI_Comm comm,DM *dm2)
153 {
154   PetscErrorCode          ierr;
155   PetscContainer          isnes;
156   DM_SNESVI               *dmsnesvi1;
157   Vec                     upper,lower,values,F;
158   IS                      inactive;
159   VecScatter              inject;
160 
161   PetscFunctionBegin;
162   ierr = PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject *)&isnes);CHKERRQ(ierr);
163   if (!isnes) SETERRQ(((PetscObject)dm1)->comm,PETSC_ERR_PLIB,"Composed VI data structure is missing");
164   ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);CHKERRQ(ierr);
165 
166   /* get the original coarsen */
167   ierr = (*dmsnesvi1->coarsen)(dm1,comm,dm2);CHKERRQ(ierr);
168 
169   /* not sure why this extra reference is needed, but without the dm2 disappears too early */
170   ierr = PetscObjectReference((PetscObject)*dm2);CHKERRQ(ierr);
171 
172   /* need to set back global vectors in order to use the original injection */
173   ierr = DMClearGlobalVectors(dm1);CHKERRQ(ierr);
174   dm1->ops->createglobalvector = dmsnesvi1->createglobalvector;
175   ierr = DMGetInjection(*dm2,dm1,&inject);CHKERRQ(ierr);
176   ierr = DMCreateGlobalVector(*dm2,&upper);CHKERRQ(ierr);
177   ierr = VecScatterBegin(inject,dmsnesvi1->upper,upper,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
178   ierr = VecScatterEnd(inject,dmsnesvi1->upper,upper,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
179   ierr = DMCreateGlobalVector(*dm2,&lower);CHKERRQ(ierr);
180   ierr = VecScatterBegin(inject,dmsnesvi1->lower,lower,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
181   ierr = VecScatterEnd(inject,dmsnesvi1->lower,lower,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
182   ierr = DMCreateGlobalVector(*dm2,&values);CHKERRQ(ierr);
183   ierr = VecScatterBegin(inject,dmsnesvi1->values,values,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
184   ierr = VecScatterEnd(inject,dmsnesvi1->values,values,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
185   ierr = DMCreateGlobalVector(*dm2,&F);CHKERRQ(ierr);
186   ierr = VecScatterBegin(inject,dmsnesvi1->F,F,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
187   ierr = VecScatterEnd(inject,dmsnesvi1->F,F,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
188   ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
189   ierr = DMClearGlobalVectors(dm1);CHKERRQ(ierr);
190   dm1->ops->createglobalvector = DMCreateGlobalVector_SNESVI;
191   ierr = SNESVIGetInActiveSetIS(upper,lower,values,F,&inactive);CHKERRQ(ierr);
192   ierr = DMSetVI(*dm2,upper,lower,values,F,inactive);CHKERRQ(ierr);
193   ierr = VecDestroy(&upper);CHKERRQ(ierr);
194   ierr = VecDestroy(&lower);CHKERRQ(ierr);
195   ierr = VecDestroy(&values);CHKERRQ(ierr);
196   ierr = VecDestroy(&F);CHKERRQ(ierr);
197   ierr = ISDestroy(&inactive);CHKERRQ(ierr);
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "DMDestroy_SNESVI"
203 PetscErrorCode DMDestroy_SNESVI(DM_SNESVI *dmsnesvi)
204 {
205   PetscErrorCode ierr;
206 
207   PetscFunctionBegin;
208   /* reset the base methods in the DM object that were changed when the DM_SNESVI was reset */
209   dmsnesvi->dm->ops->getinterpolation   = dmsnesvi->getinterpolation;
210   dmsnesvi->dm->ops->coarsen            = dmsnesvi->coarsen;
211   dmsnesvi->dm->ops->createglobalvector = dmsnesvi->createglobalvector;
212   /* need to clear out this vectors because some of them may not have a reference to the DM
213     but they are counted as having references to the DM in DMDestroy() */
214   ierr = DMClearGlobalVectors(dmsnesvi->dm);CHKERRQ(ierr);
215 
216   ierr = VecDestroy(&dmsnesvi->upper);CHKERRQ(ierr);
217   ierr = VecDestroy(&dmsnesvi->lower);CHKERRQ(ierr);
218   ierr = VecDestroy(&dmsnesvi->values);CHKERRQ(ierr);
219   ierr = VecDestroy(&dmsnesvi->F);CHKERRQ(ierr);
220   ierr = ISDestroy(&dmsnesvi->inactive);CHKERRQ(ierr);
221   ierr = PetscFree(dmsnesvi);CHKERRQ(ierr);
222   PetscFunctionReturn(0);
223 }
224 
225 #undef __FUNCT__
226 #define __FUNCT__ "DMSetVI"
227 /*
228      DMSetVI - Marks a DM as associated with a VI problem. This causes the interpolation/restriction operators to
229                be restricted to only those variables NOT associated with active constraints.
230 
231 */
232 PetscErrorCode  DMSetVI(DM dm,Vec upper,Vec lower,Vec values,Vec F,IS inactive)
233 {
234   PetscErrorCode          ierr;
235   PetscContainer          isnes;
236   DM_SNESVI               *dmsnesvi;
237 
238   PetscFunctionBegin;
239   if (!dm) PetscFunctionReturn(0);
240 
241   ierr = PetscObjectReference((PetscObject)upper);CHKERRQ(ierr);
242   ierr = PetscObjectReference((PetscObject)lower);CHKERRQ(ierr);
243   ierr = PetscObjectReference((PetscObject)values);CHKERRQ(ierr);
244   ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
245   ierr = PetscObjectReference((PetscObject)inactive);CHKERRQ(ierr);
246 
247   ierr = PetscObjectQuery((PetscObject)dm,"VI",(PetscObject *)&isnes);CHKERRQ(ierr);
248   if (!isnes) {
249     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&isnes);CHKERRQ(ierr);
250     ierr = PetscContainerSetUserDestroy(isnes,(PetscErrorCode (*)(void*))DMDestroy_SNESVI);CHKERRQ(ierr);
251     ierr = PetscNew(DM_SNESVI,&dmsnesvi);CHKERRQ(ierr);
252     ierr = PetscContainerSetPointer(isnes,(void*)dmsnesvi);CHKERRQ(ierr);
253     ierr = PetscObjectCompose((PetscObject)dm,"VI",(PetscObject)isnes);CHKERRQ(ierr);
254     ierr = PetscContainerDestroy(&isnes);CHKERRQ(ierr);
255     dmsnesvi->getinterpolation   = dm->ops->getinterpolation;
256     dm->ops->getinterpolation    = DMGetInterpolation_SNESVI;
257     dmsnesvi->coarsen            = dm->ops->coarsen;
258     dm->ops->coarsen             = DMCoarsen_SNESVI;
259     dmsnesvi->createglobalvector = dm->ops->createglobalvector;
260     dm->ops->createglobalvector  = DMCreateGlobalVector_SNESVI;
261     /* since these vectors may reference the DM, need to remove circle referencing */
262     ierr = PetscObjectRemoveReference((PetscObject)upper,"DM");CHKERRQ(ierr);
263     ierr = PetscObjectRemoveReference((PetscObject)lower,"DM");CHKERRQ(ierr);
264     ierr = PetscObjectRemoveReference((PetscObject)values,"DM");CHKERRQ(ierr);
265     ierr = PetscObjectRemoveReference((PetscObject)F,"DM");CHKERRQ(ierr);
266   } else {
267     ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi);CHKERRQ(ierr);
268     ierr = VecDestroy(&dmsnesvi->upper);CHKERRQ(ierr);
269     ierr = VecDestroy(&dmsnesvi->lower);CHKERRQ(ierr);
270     ierr = VecDestroy(&dmsnesvi->values);CHKERRQ(ierr);
271     ierr = VecDestroy(&dmsnesvi->F);CHKERRQ(ierr);
272     ierr = ISDestroy(&dmsnesvi->inactive);CHKERRQ(ierr);
273   }
274   ierr = DMClearGlobalVectors(dm);CHKERRQ(ierr);
275   ierr = ISGetLocalSize(inactive,&dmsnesvi->n);CHKERRQ(ierr);
276   dmsnesvi->upper    = upper;
277   dmsnesvi->lower    = lower;
278   dmsnesvi->values   = values;
279   dmsnesvi->F        = F;
280   dmsnesvi->inactive = inactive;
281   dmsnesvi->dm       = dm;
282   PetscFunctionReturn(0);
283 }
284 
285 #undef __FUNCT__
286 #define __FUNCT__ "DMDestroyVI"
287 /*
288      DMDestroyVI - Frees the DM_SNESVI object contained in the DM
289          - also resets the function pointers in the DM for getinterpolation() etc to use the original DM
290 */
291 PetscErrorCode  DMDestroyVI(DM dm)
292 {
293   PetscErrorCode          ierr;
294 
295   PetscFunctionBegin;
296   if (!dm) PetscFunctionReturn(0);
297   ierr = PetscObjectCompose((PetscObject)dm,"VI",(PetscObject)PETSC_NULL);CHKERRQ(ierr);
298   PetscFunctionReturn(0);
299 }
300 
301 /* --------------------------------------------------------------------------------------------------------*/
302 
303 #undef __FUNCT__
304 #define __FUNCT__ "SNESMonitorVI"
305 PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
306 {
307   PetscErrorCode    ierr;
308   SNES_VI            *vi = (SNES_VI*)snes->data;
309   PetscViewer        viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
310   const PetscScalar  *x,*xl,*xu,*f;
311   PetscInt           i,n,act = 0;
312   PetscReal          rnorm,fnorm;
313 
314   PetscFunctionBegin;
315   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
316   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
317   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
318   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
319   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
320 
321   rnorm = 0.0;
322   for (i=0; i<n; i++) {
323     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
324     else act++;
325   }
326   ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
327   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
328   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
329   ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
330   ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
331   fnorm = sqrt(fnorm);
332   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
333   ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES VI Function norm %14.12e Active constraints %D\n",its,fnorm,act);CHKERRQ(ierr);
334   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
335   PetscFunctionReturn(0);
336 }
337 
338 /*
339      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
340     || 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
341     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
342     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
343 */
344 #undef __FUNCT__
345 #define __FUNCT__ "SNESVICheckLocalMin_Private"
346 PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
347 {
348   PetscReal      a1;
349   PetscErrorCode ierr;
350   PetscBool     hastranspose;
351 
352   PetscFunctionBegin;
353   *ismin = PETSC_FALSE;
354   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
355   if (hastranspose) {
356     /* Compute || J^T F|| */
357     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
358     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
359     ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr);
360     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
361   } else {
362     Vec         work;
363     PetscScalar result;
364     PetscReal   wnorm;
365 
366     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
367     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
368     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
369     ierr = MatMult(A,W,work);CHKERRQ(ierr);
370     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
371     ierr = VecDestroy(&work);CHKERRQ(ierr);
372     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
373     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr);
374     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
375   }
376   PetscFunctionReturn(0);
377 }
378 
379 /*
380      Checks if J^T(F - J*X) = 0
381 */
382 #undef __FUNCT__
383 #define __FUNCT__ "SNESVICheckResidual_Private"
384 PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
385 {
386   PetscReal      a1,a2;
387   PetscErrorCode ierr;
388   PetscBool     hastranspose;
389 
390   PetscFunctionBegin;
391   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
392   if (hastranspose) {
393     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
394     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
395 
396     /* Compute || J^T W|| */
397     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
398     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
399     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
400     if (a1 != 0.0) {
401       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr);
402     }
403   }
404   PetscFunctionReturn(0);
405 }
406 
407 /*
408   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
409 
410   Notes:
411   The convergence criterion currently implemented is
412   merit < abstol
413   merit < rtol*merit_initial
414 */
415 #undef __FUNCT__
416 #define __FUNCT__ "SNESDefaultConverged_VI"
417 PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
418 {
419   PetscErrorCode ierr;
420 
421   PetscFunctionBegin;
422   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
423   PetscValidPointer(reason,6);
424 
425   *reason = SNES_CONVERGED_ITERATING;
426 
427   if (!it) {
428     /* set parameter for default relative tolerance convergence test */
429     snes->ttol = fnorm*snes->rtol;
430   }
431   if (fnorm != fnorm) {
432     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
433     *reason = SNES_DIVERGED_FNORM_NAN;
434   } else if (fnorm < snes->abstol) {
435     ierr = PetscInfo2(snes,"Converged due to function norm %G < %G\n",fnorm,snes->abstol);CHKERRQ(ierr);
436     *reason = SNES_CONVERGED_FNORM_ABS;
437   } else if (snes->nfuncs >= snes->max_funcs) {
438     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
439     *reason = SNES_DIVERGED_FUNCTION_COUNT;
440   }
441 
442   if (it && !*reason) {
443     if (fnorm < snes->ttol) {
444       ierr = PetscInfo2(snes,"Converged due to function norm %G < %G (relative tolerance)\n",fnorm,snes->ttol);CHKERRQ(ierr);
445       *reason = SNES_CONVERGED_FNORM_RELATIVE;
446     }
447   }
448   PetscFunctionReturn(0);
449 }
450 
451 /*
452   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
453 
454   Input Parameter:
455 . phi - the semismooth function
456 
457   Output Parameter:
458 . merit - the merit function
459 . phinorm - ||phi||
460 
461   Notes:
462   The merit function for the mixed complementarity problem is defined as
463      merit = 0.5*phi^T*phi
464 */
465 #undef __FUNCT__
466 #define __FUNCT__ "SNESVIComputeMeritFunction"
467 static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm)
468 {
469   PetscErrorCode ierr;
470 
471   PetscFunctionBegin;
472   ierr = VecNormBegin(phi,NORM_2,phinorm);CHKERRQ(ierr);
473   ierr = VecNormEnd(phi,NORM_2,phinorm);CHKERRQ(ierr);
474 
475   *merit = 0.5*(*phinorm)*(*phinorm);
476   PetscFunctionReturn(0);
477 }
478 
479 PETSC_STATIC_INLINE PetscScalar Phi(PetscScalar a,PetscScalar b)
480 {
481   return a + b - sqrt(a*a + b*b);
482 }
483 
484 PETSC_STATIC_INLINE PetscScalar DPhi(PetscScalar a,PetscScalar b)
485 {
486   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return  1.0 - a/ sqrt(a*a + b*b);
487   else return .5;
488 }
489 
490 /*
491    SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form.
492 
493    Input Parameters:
494 .  snes - the SNES context
495 .  x - current iterate
496 .  functx - user defined function context
497 
498    Output Parameters:
499 .  phi - Semismooth function
500 
501 */
502 #undef __FUNCT__
503 #define __FUNCT__ "SNESVIComputeFunction"
504 static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx)
505 {
506   PetscErrorCode  ierr;
507   SNES_VI         *vi = (SNES_VI*)snes->data;
508   Vec             Xl = vi->xl,Xu = vi->xu,F = snes->vec_func;
509   PetscScalar     *phi_arr,*x_arr,*f_arr,*l,*u;
510   PetscInt        i,nlocal;
511 
512   PetscFunctionBegin;
513   ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr);
514   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
515   ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr);
516   ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr);
517   ierr = VecGetArray(Xl,&l);CHKERRQ(ierr);
518   ierr = VecGetArray(Xu,&u);CHKERRQ(ierr);
519   ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr);
520 
521   for (i=0;i < nlocal;i++) {
522     if ((PetscRealPart(l[i]) <= SNES_VI_NINF) && (PetscRealPart(u[i]) >= SNES_VI_INF)) { /* no constraints on variable */
523       phi_arr[i] = f_arr[i];
524     } else if (PetscRealPart(l[i]) <= SNES_VI_NINF) {                      /* upper bound on variable only */
525       phi_arr[i] = -Phi(u[i] - x_arr[i],-f_arr[i]);
526     } else if (PetscRealPart(u[i]) >= SNES_VI_INF) {                       /* lower bound on variable only */
527       phi_arr[i] = Phi(x_arr[i] - l[i],f_arr[i]);
528     } else if (l[i] == u[i]) {
529       phi_arr[i] = l[i] - x_arr[i];
530     } else {                                                /* both bounds on variable */
531       phi_arr[i] = Phi(x_arr[i] - l[i],-Phi(u[i] - x_arr[i],-f_arr[i]));
532     }
533   }
534 
535   ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr);
536   ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr);
537   ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr);
538   ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr);
539   ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr);
540   PetscFunctionReturn(0);
541 }
542 
543 /*
544    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
545                                           the semismooth jacobian.
546 */
547 #undef __FUNCT__
548 #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors"
549 PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db)
550 {
551   PetscErrorCode ierr;
552   SNES_VI      *vi = (SNES_VI*)snes->data;
553   PetscScalar    *l,*u,*x,*f,*da,*db,da1,da2,db1,db2;
554   PetscInt       i,nlocal;
555 
556   PetscFunctionBegin;
557 
558   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
559   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
560   ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr);
561   ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr);
562   ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
563   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
564   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
565 
566   for (i=0;i< nlocal;i++) {
567     if ((PetscRealPart(l[i]) <= SNES_VI_NINF) && (PetscRealPart(u[i]) >= SNES_VI_INF)) {/* no constraints on variable */
568       da[i] = 0;
569       db[i] = 1;
570     } else if (PetscRealPart(l[i]) <= SNES_VI_NINF) {                     /* upper bound on variable only */
571       da[i] = DPhi(u[i] - x[i], -f[i]);
572       db[i] = DPhi(-f[i],u[i] - x[i]);
573     } else if (PetscRealPart(u[i]) >= SNES_VI_INF) {                      /* lower bound on variable only */
574       da[i] = DPhi(x[i] - l[i], f[i]);
575       db[i] = DPhi(f[i],x[i] - l[i]);
576     } else if (l[i] == u[i]) {                              /* fixed variable */
577       da[i] = 1;
578       db[i] = 0;
579     } else {                                                /* upper and lower bounds on variable */
580       da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
581       db1 = DPhi(-Phi(u[i] - x[i], -f[i]),x[i] - l[i]);
582       da2 = DPhi(u[i] - x[i], -f[i]);
583       db2 = DPhi(-f[i],u[i] - x[i]);
584       da[i] = da1 + db1*da2;
585       db[i] = db1*db2;
586     }
587   }
588 
589   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
590   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
591   ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr);
592   ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr);
593   ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
594   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
595   PetscFunctionReturn(0);
596 }
597 
598 /*
599    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.
600 
601    Input Parameters:
602 .  Da       - Diagonal shift vector for the semismooth jacobian.
603 .  Db       - Row scaling vector for the semismooth jacobian.
604 
605    Output Parameters:
606 .  jac      - semismooth jacobian
607 .  jac_pre  - optional preconditioning matrix
608 
609    Notes:
610    The semismooth jacobian matrix is given by
611    jac = Da + Db*jacfun
612    where Db is the row scaling matrix stored as a vector,
613          Da is the diagonal perturbation matrix stored as a vector
614    and   jacfun is the jacobian of the original nonlinear function.
615 */
616 #undef __FUNCT__
617 #define __FUNCT__ "SNESVIComputeJacobian"
618 PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db)
619 {
620   PetscErrorCode ierr;
621 
622   /* Do row scaling  and add diagonal perturbation */
623   ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr);
624   ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr);
625   if (jac != jac_pre) { /* If jac and jac_pre are different */
626     ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL);
627     ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr);
628   }
629   PetscFunctionReturn(0);
630 }
631 
632 /*
633    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
634 
635    Input Parameters:
636    phi - semismooth function.
637    H   - semismooth jacobian
638 
639    Output Parameters:
640    dpsi - merit function gradient
641 
642    Notes:
643   The merit function gradient is computed as follows
644         dpsi = H^T*phi
645 */
646 #undef __FUNCT__
647 #define __FUNCT__ "SNESVIComputeMeritFunctionGradient"
648 PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi)
649 {
650   PetscErrorCode ierr;
651 
652   PetscFunctionBegin;
653   ierr = MatMultTranspose(H,phi,dpsi);CHKERRQ(ierr);
654   PetscFunctionReturn(0);
655 }
656 
657 /* -------------------------------------------------------------------------- */
658 /*
659    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
660 
661    Input Parameters:
662 .  SNES - nonlinear solver context
663 
664    Output Parameters:
665 .  X - Bound projected X
666 
667 */
668 
669 #undef __FUNCT__
670 #define __FUNCT__ "SNESVIProjectOntoBounds"
671 PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
672 {
673   PetscErrorCode    ierr;
674   SNES_VI           *vi = (SNES_VI*)snes->data;
675   const PetscScalar *xl,*xu;
676   PetscScalar       *x;
677   PetscInt          i,n;
678 
679   PetscFunctionBegin;
680   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
681   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
682   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
683   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
684 
685   for(i = 0;i<n;i++) {
686     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
687     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
688   }
689   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
690   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
691   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
692   PetscFunctionReturn(0);
693 }
694 
695 /*  --------------------------------------------------------------------
696 
697      This file implements a semismooth truncated Newton method with a line search,
698      for solving a system of nonlinear equations in complementarity form, using the KSP, Vec,
699      and Mat interfaces for linear solvers, vectors, and matrices,
700      respectively.
701 
702      The following basic routines are required for each nonlinear solver:
703           SNESCreate_XXX()          - Creates a nonlinear solver context
704           SNESSetFromOptions_XXX()  - Sets runtime options
705           SNESSolve_XXX()           - Solves the nonlinear system
706           SNESDestroy_XXX()         - Destroys the nonlinear solver context
707      The suffix "_XXX" denotes a particular implementation, in this case
708      we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving
709      systems of nonlinear equations with a line search (LS) method.
710      These routines are actually called via the common user interface
711      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
712      SNESDestroy(), so the application code interface remains identical
713      for all nonlinear solvers.
714 
715      Another key routine is:
716           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
717      by setting data structures and options.   The interface routine SNESSetUp()
718      is not usually called directly by the user, but instead is called by
719      SNESSolve() if necessary.
720 
721      Additional basic routines are:
722           SNESView_XXX()            - Prints details of runtime options that
723                                       have actually been used.
724      These are called by application codes via the interface routines
725      SNESView().
726 
727      The various types of solvers (preconditioners, Krylov subspace methods,
728      nonlinear solvers, timesteppers) are all organized similarly, so the
729      above description applies to these categories also.
730 
731     -------------------------------------------------------------------- */
732 /*
733    SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton
734    method using a line search.
735 
736    Input Parameters:
737 .  snes - the SNES context
738 
739    Output Parameter:
740 .  outits - number of iterations until termination
741 
742    Application Interface Routine: SNESSolve()
743 
744    Notes:
745    This implements essentially a semismooth Newton method with a
746    line search. The default line search does not do any line seach
747    but rather takes a full newton step.
748 */
749 #undef __FUNCT__
750 #define __FUNCT__ "SNESSolveVI_SS"
751 PetscErrorCode SNESSolveVI_SS(SNES snes)
752 {
753   SNES_VI            *vi = (SNES_VI*)snes->data;
754   PetscErrorCode     ierr;
755   PetscInt           maxits,i,lits;
756   PetscBool          lssucceed;
757   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
758   PetscReal          gnorm,xnorm=0,ynorm;
759   Vec                Y,X,F,G,W;
760   KSPConvergedReason kspreason;
761 
762   PetscFunctionBegin;
763   vi->computeuserfunction    = snes->ops->computefunction;
764   snes->ops->computefunction = SNESVIComputeFunction;
765 
766   snes->numFailures            = 0;
767   snes->numLinearSolveFailures = 0;
768   snes->reason                 = SNES_CONVERGED_ITERATING;
769 
770   maxits	= snes->max_its;	/* maximum number of iterations */
771   X		= snes->vec_sol;	/* solution vector */
772   F		= snes->vec_func;	/* residual vector */
773   Y		= snes->work[0];	/* work vectors */
774   G		= snes->work[1];
775   W		= snes->work[2];
776 
777   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
778   snes->iter = 0;
779   snes->norm = 0.0;
780   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
781 
782   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
783   ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr);
784   if (snes->domainerror) {
785     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
786     snes->ops->computefunction = vi->computeuserfunction;
787     PetscFunctionReturn(0);
788   }
789    /* Compute Merit function */
790   ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr);
791 
792   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
793   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
794   if (PetscIsInfOrNanReal(vi->merit)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
795 
796   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
797   snes->norm = vi->phinorm;
798   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
799   SNESLogConvHistory(snes,vi->phinorm,0);
800   ierr = SNESMonitor(snes,0,vi->phinorm);CHKERRQ(ierr);
801 
802   /* set parameter for default relative tolerance convergence test */
803   snes->ttol = vi->phinorm*snes->rtol;
804   /* test convergence */
805   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
806   if (snes->reason) {
807     snes->ops->computefunction = vi->computeuserfunction;
808     PetscFunctionReturn(0);
809   }
810 
811   for (i=0; i<maxits; i++) {
812 
813     /* Call general purpose update function */
814     if (snes->ops->update) {
815       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
816     }
817 
818     /* Solve J Y = Phi, where J is the semismooth jacobian */
819     /* Get the nonlinear function jacobian */
820     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
821     /* Get the diagonal shift and row scaling vectors */
822     ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr);
823     /* Compute the semismooth jacobian */
824     ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr);
825     /* Compute the merit function gradient */
826     ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr);
827     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
828     ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr);
829     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
830 
831     if (kspreason < 0) {
832       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
833         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
834         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
835         break;
836       }
837     }
838     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
839     snes->linear_its += lits;
840     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
841     /*
842     if (vi->precheckstep) {
843       PetscBool changed_y = PETSC_FALSE;
844       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
845     }
846 
847     if (PetscLogPrintInfo){
848       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
849     }
850     */
851     /* Compute a (scaled) negative update in the line search routine:
852          Y <- X - lambda*Y
853        and evaluate G = function(Y) (depends on the line search).
854     */
855     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
856     ynorm = 1; gnorm = vi->phinorm;
857     ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
858     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
859     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
860     if (snes->domainerror) {
861       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
862       snes->ops->computefunction = vi->computeuserfunction;
863       PetscFunctionReturn(0);
864     }
865     if (!lssucceed) {
866       if (++snes->numFailures >= snes->maxFailures) {
867 	PetscBool ismin;
868         snes->reason = SNES_DIVERGED_LINE_SEARCH;
869         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
870         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
871         break;
872       }
873     }
874     /* Update function and solution vectors */
875     vi->phinorm = gnorm;
876     vi->merit = 0.5*vi->phinorm*vi->phinorm;
877     ierr = VecCopy(G,vi->phi);CHKERRQ(ierr);
878     ierr = VecCopy(W,X);CHKERRQ(ierr);
879     /* Monitor convergence */
880     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
881     snes->iter = i+1;
882     snes->norm = vi->phinorm;
883     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
884     SNESLogConvHistory(snes,snes->norm,lits);
885     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
886     /* Test for convergence, xnorm = || X || */
887     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
888     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
889     if (snes->reason) break;
890   }
891   if (i == maxits) {
892     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
893     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
894   }
895   snes->ops->computefunction = vi->computeuserfunction;
896   PetscFunctionReturn(0);
897 }
898 
899 #undef __FUNCT__
900 #define __FUNCT__ "SNESVIGetActiveSetIS"
901 /*
902    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
903 
904    Input parameter
905 .  snes - the SNES context
906 .  X    - the snes solution vector
907 .  F    - the nonlinear function vector
908 
909    Output parameter
910 .  ISact - active set index set
911  */
912 PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS* ISact)
913 {
914   PetscErrorCode   ierr;
915   SNES_VI          *vi = (SNES_VI*)snes->data;
916   Vec               Xl=vi->xl,Xu=vi->xu;
917   const PetscScalar *x,*f,*xl,*xu;
918   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
919 
920   PetscFunctionBegin;
921   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
922   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
923   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
924   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
925   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
926   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
927   /* Compute active set size */
928   for (i=0; i < nlocal;i++) {
929     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) nloc_isact++;
930   }
931 
932   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
933 
934   /* Set active set indices */
935   for(i=0; i < nlocal; i++) {
936     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) idx_act[i1++] = ilow+i;
937   }
938 
939    /* Create active set IS */
940   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
941 
942   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
943   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
944   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
945   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
946   PetscFunctionReturn(0);
947 }
948 
949 #undef __FUNCT__
950 #define __FUNCT__ "SNESVICreateIndexSets_RS"
951 PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS* ISact,IS* ISinact)
952 {
953   PetscErrorCode     ierr;
954 
955   PetscFunctionBegin;
956   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
957   ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr);
958   PetscFunctionReturn(0);
959 }
960 
961 /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */
962 #undef __FUNCT__
963 #define __FUNCT__ "SNESVICreateSubVectors"
964 PetscErrorCode SNESVICreateSubVectors(SNES snes,PetscInt n,Vec* newv)
965 {
966   PetscErrorCode ierr;
967   Vec            v;
968 
969   PetscFunctionBegin;
970   ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr);
971   ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr);
972   ierr = VecSetFromOptions(v);CHKERRQ(ierr);
973   *newv = v;
974 
975   PetscFunctionReturn(0);
976 }
977 
978 /* Resets the snes PC and KSP when the active set sizes change */
979 #undef __FUNCT__
980 #define __FUNCT__ "SNESVIResetPCandKSP"
981 PetscErrorCode SNESVIResetPCandKSP(SNES snes,Mat Amat,Mat Pmat)
982 {
983   PetscErrorCode         ierr;
984   KSP                    snesksp;
985 
986   PetscFunctionBegin;
987   ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr);
988   ierr = KSPReset(snesksp);CHKERRQ(ierr);
989 
990   /*
991   KSP                    kspnew;
992   PC                     pcnew;
993   const MatSolverPackage stype;
994 
995 
996   ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr);
997   kspnew->pc_side = snesksp->pc_side;
998   kspnew->rtol    = snesksp->rtol;
999   kspnew->abstol    = snesksp->abstol;
1000   kspnew->max_it  = snesksp->max_it;
1001   ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr);
1002   ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr);
1003   ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr);
1004   ierr = PCSetOperators(kspnew->pc,Amat,Pmat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
1005   ierr = PCFactorGetMatSolverPackage(snesksp->pc,&stype);CHKERRQ(ierr);
1006   ierr = PCFactorSetMatSolverPackage(kspnew->pc,stype);CHKERRQ(ierr);
1007   ierr = KSPDestroy(&snesksp);CHKERRQ(ierr);
1008   snes->ksp = kspnew;
1009   ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr);
1010    ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr);*/
1011   PetscFunctionReturn(0);
1012 }
1013 
1014 
1015 #undef __FUNCT__
1016 #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
1017 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X,PetscReal *fnorm)
1018 {
1019   PetscErrorCode    ierr;
1020   SNES_VI           *vi = (SNES_VI*)snes->data;
1021   const PetscScalar *x,*xl,*xu,*f;
1022   PetscInt          i,n;
1023   PetscReal         rnorm;
1024 
1025   PetscFunctionBegin;
1026   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
1027   ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr);
1028   ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr);
1029   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1030   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
1031   rnorm = 0.0;
1032   for (i=0; i<n; i++) {
1033     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
1034   }
1035   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
1036   ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr);
1037   ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr);
1038   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
1039   ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
1040   *fnorm = sqrt(*fnorm);
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 /* Variational Inequality solver using reduce space method. No semismooth algorithm is
1045    implemented in this algorithm. It basically identifies the active constraints and does
1046    a linear solve on the other variables (those not associated with the active constraints). */
1047 #undef __FUNCT__
1048 #define __FUNCT__ "SNESSolveVI_RS"
1049 PetscErrorCode SNESSolveVI_RS(SNES snes)
1050 {
1051   SNES_VI          *vi = (SNES_VI*)snes->data;
1052   PetscErrorCode    ierr;
1053   PetscInt          maxits,i,lits;
1054   PetscBool         lssucceed;
1055   MatStructure      flg = DIFFERENT_NONZERO_PATTERN;
1056   PetscReal         fnorm,gnorm,xnorm=0,ynorm;
1057   Vec                Y,X,F,G,W;
1058   KSPConvergedReason kspreason;
1059 
1060   PetscFunctionBegin;
1061   snes->numFailures            = 0;
1062   snes->numLinearSolveFailures = 0;
1063   snes->reason                 = SNES_CONVERGED_ITERATING;
1064 
1065   maxits	= snes->max_its;	/* maximum number of iterations */
1066   X		= snes->vec_sol;	/* solution vector */
1067   F		= snes->vec_func;	/* residual vector */
1068   Y		= snes->work[0];	/* work vectors */
1069   G		= snes->work[1];
1070   W		= snes->work[2];
1071 
1072   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1073   snes->iter = 0;
1074   snes->norm = 0.0;
1075   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1076 
1077   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
1078   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1079   if (snes->domainerror) {
1080     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1081     PetscFunctionReturn(0);
1082   }
1083   ierr = SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);CHKERRQ(ierr);
1084   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
1085   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
1086   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1087 
1088   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1089   snes->norm = fnorm;
1090   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1091   SNESLogConvHistory(snes,fnorm,0);
1092   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1093 
1094   /* set parameter for default relative tolerance convergence test */
1095   snes->ttol = fnorm*snes->rtol;
1096   /* test convergence */
1097   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1098   if (snes->reason) PetscFunctionReturn(0);
1099 
1100   for (i=0; i<maxits; i++) {
1101 
1102     IS         IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
1103     IS         IS_redact; /* redundant active set */
1104     VecScatter scat_act,scat_inact;
1105     PetscInt   nis_act,nis_inact;
1106     Vec        Y_act,Y_inact,F_inact;
1107     Mat        jac_inact_inact,prejac_inact_inact;
1108     IS         keptrows;
1109     PetscBool  isequal;
1110 
1111     /* Call general purpose update function */
1112     if (snes->ops->update) {
1113       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1114     }
1115     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
1116 
1117     /* Create active and inactive index sets */
1118 
1119     /*original
1120     ierr = SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&IS_inact);CHKERRQ(ierr);
1121      */
1122     ierr = SNESVIGetActiveSetIS(snes,X,F,&IS_act);CHKERRQ(ierr);
1123 
1124     if (vi->checkredundancy) {
1125       (*vi->checkredundancy)(snes,IS_act,&IS_redact,vi->ctxP);CHKERRQ(ierr);
1126       ierr = ISComplement(IS_redact,X->map->rstart,X->map->rend,&IS_inact);CHKERRQ(ierr);
1127       ierr = ISDestroy(&IS_redact);CHKERRQ(ierr);
1128     } else {
1129       ierr = ISComplement(IS_act,X->map->rstart,X->map->rend,&IS_inact);CHKERRQ(ierr);
1130     }
1131 
1132 
1133     /* Create inactive set submatrix */
1134     ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
1135 
1136     ierr = MatFindNonzeroRows(jac_inact_inact,&keptrows);CHKERRQ(ierr);
1137     if (keptrows) {
1138       PetscInt       cnt,*nrows,k;
1139       const PetscInt *krows,*inact;
1140       PetscInt       rstart=jac_inact_inact->rmap->rstart;
1141 
1142       ierr = MatDestroy(&jac_inact_inact);CHKERRQ(ierr);
1143       ierr = ISDestroy(&IS_act);CHKERRQ(ierr);
1144 
1145       ierr = ISGetLocalSize(keptrows,&cnt);CHKERRQ(ierr);
1146       ierr = ISGetIndices(keptrows,&krows);CHKERRQ(ierr);
1147       ierr = ISGetIndices(IS_inact,&inact);CHKERRQ(ierr);
1148       ierr = PetscMalloc(cnt*sizeof(PetscInt),&nrows);CHKERRQ(ierr);
1149       for (k=0; k<cnt; k++) {
1150         nrows[k] = inact[krows[k]-rstart];
1151       }
1152       ierr = ISRestoreIndices(keptrows,&krows);CHKERRQ(ierr);
1153       ierr = ISRestoreIndices(IS_inact,&inact);CHKERRQ(ierr);
1154       ierr = ISDestroy(&keptrows);CHKERRQ(ierr);
1155       ierr = ISDestroy(&IS_inact);CHKERRQ(ierr);
1156 
1157       ierr = ISCreateGeneral(PETSC_COMM_WORLD,cnt,nrows,PETSC_OWN_POINTER,&IS_inact);CHKERRQ(ierr);
1158       ierr = ISComplement(IS_inact,F->map->rstart,F->map->rend,&IS_act);CHKERRQ(ierr);
1159       ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
1160     }
1161     ierr = DMSetVI(snes->dm,vi->xu,vi->xl,X,F,IS_inact);CHKERRQ(ierr);
1162 
1163     /* Get sizes of active and inactive sets */
1164     ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr);
1165     ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr);
1166 
1167     /* Create active and inactive set vectors */
1168     ierr = SNESVICreateSubVectors(snes,nis_inact,&F_inact);CHKERRQ(ierr);
1169     ierr = SNESVICreateSubVectors(snes,nis_act,&Y_act);CHKERRQ(ierr);
1170     ierr = SNESVICreateSubVectors(snes,nis_inact,&Y_inact);CHKERRQ(ierr);
1171 
1172     /* Create scatter contexts */
1173     ierr = VecScatterCreate(Y,IS_act,Y_act,PETSC_NULL,&scat_act);CHKERRQ(ierr);
1174     ierr = VecScatterCreate(Y,IS_inact,Y_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr);
1175 
1176     /* Do a vec scatter to active and inactive set vectors */
1177     ierr = VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1178     ierr = VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1179 
1180     ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1181     ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1182 
1183     ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1184     ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1185 
1186     /* Active set direction = 0 */
1187     ierr = VecSet(Y_act,0);CHKERRQ(ierr);
1188     if (snes->jacobian != snes->jacobian_pre) {
1189       ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr);
1190     } else prejac_inact_inact = jac_inact_inact;
1191 
1192     ierr = ISEqual(vi->IS_inact_prev,IS_inact,&isequal);CHKERRQ(ierr);
1193     if (!isequal) {
1194       ierr = SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);CHKERRQ(ierr);
1195       flg  = DIFFERENT_NONZERO_PATTERN;
1196     }
1197 
1198     /*      ierr = ISView(IS_inact,0);CHKERRQ(ierr); */
1199     /*      ierr = ISView(IS_act,0);CHKERRQ(ierr);*/
1200     /*      ierr = MatView(snes->jacobian_pre,0); */
1201 
1202     ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr);
1203     ierr = KSPSetUp(snes->ksp);CHKERRQ(ierr);
1204     {
1205       PC        pc;
1206       PetscBool flg;
1207       ierr = KSPGetPC(snes->ksp,&pc);CHKERRQ(ierr);
1208       ierr = PetscTypeCompare((PetscObject)pc,PCFIELDSPLIT,&flg);CHKERRQ(ierr);
1209       if (flg) {
1210         KSP      *subksps;
1211         ierr = PCFieldSplitGetSubKSP(pc,PETSC_NULL,&subksps);CHKERRQ(ierr);
1212         ierr = KSPGetPC(subksps[0],&pc);CHKERRQ(ierr);
1213         ierr = PetscFree(subksps);CHKERRQ(ierr);
1214         ierr = PetscTypeCompare((PetscObject)pc,PCBJACOBI,&flg);CHKERRQ(ierr);
1215         if (flg) {
1216           PetscInt       n,N = 101*101,j,cnts[3] = {0,0,0};
1217           const PetscInt *ii;
1218 
1219           ierr = ISGetSize(IS_inact,&n);CHKERRQ(ierr);
1220           ierr = ISGetIndices(IS_inact,&ii);CHKERRQ(ierr);
1221           for (j=0; j<n; j++) {
1222             if (ii[j] < N) cnts[0]++;
1223             else if (ii[j] < 2*N) cnts[1]++;
1224             else if (ii[j] < 3*N) cnts[2]++;
1225           }
1226           ierr = ISRestoreIndices(IS_inact,&ii);CHKERRQ(ierr);
1227 
1228           ierr = PCBJacobiSetTotalBlocks(pc,3,cnts);CHKERRQ(ierr);
1229         }
1230       }
1231     }
1232 
1233     ierr = SNES_KSPSolve(snes,snes->ksp,F_inact,Y_inact);CHKERRQ(ierr);
1234     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
1235     if (kspreason < 0) {
1236       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
1237         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
1238         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
1239         break;
1240       }
1241      }
1242 
1243     ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1244     ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1245     ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1246     ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1247 
1248     ierr = VecDestroy(&F_inact);CHKERRQ(ierr);
1249     ierr = VecDestroy(&Y_act);CHKERRQ(ierr);
1250     ierr = VecDestroy(&Y_inact);CHKERRQ(ierr);
1251     ierr = VecScatterDestroy(&scat_act);CHKERRQ(ierr);
1252     ierr = VecScatterDestroy(&scat_inact);CHKERRQ(ierr);
1253     ierr = ISDestroy(&IS_act);CHKERRQ(ierr);
1254     if (!isequal) {
1255       ierr = ISDestroy(&vi->IS_inact_prev);CHKERRQ(ierr);
1256       ierr = ISDuplicate(IS_inact,&vi->IS_inact_prev);CHKERRQ(ierr);
1257     }
1258     ierr = ISDestroy(&IS_inact);CHKERRQ(ierr);
1259     ierr = MatDestroy(&jac_inact_inact);CHKERRQ(ierr);
1260     if (snes->jacobian != snes->jacobian_pre) {
1261       ierr = MatDestroy(&prejac_inact_inact);CHKERRQ(ierr);
1262     }
1263     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
1264     snes->linear_its += lits;
1265     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
1266     /*
1267     if (vi->precheckstep) {
1268       PetscBool changed_y = PETSC_FALSE;
1269       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
1270     }
1271 
1272     if (PetscLogPrintInfo){
1273       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
1274     }
1275     */
1276     /* Compute a (scaled) negative update in the line search routine:
1277          Y <- X - lambda*Y
1278        and evaluate G = function(Y) (depends on the line search).
1279     */
1280     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
1281     ynorm = 1; gnorm = fnorm;
1282     ierr = (*vi->LineSearch)(snes,vi->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1283     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
1284     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
1285     if (snes->domainerror) {
1286       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1287       ierr = DMDestroyVI(snes->dm);CHKERRQ(ierr);
1288       PetscFunctionReturn(0);
1289     }
1290     if (!lssucceed) {
1291       if (++snes->numFailures >= snes->maxFailures) {
1292 	PetscBool ismin;
1293         snes->reason = SNES_DIVERGED_LINE_SEARCH;
1294         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
1295         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
1296         break;
1297       }
1298     }
1299     /* Update function and solution vectors */
1300     fnorm = gnorm;
1301     ierr = VecCopy(G,F);CHKERRQ(ierr);
1302     ierr = VecCopy(W,X);CHKERRQ(ierr);
1303     /* Monitor convergence */
1304     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1305     snes->iter = i+1;
1306     snes->norm = fnorm;
1307     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1308     SNESLogConvHistory(snes,snes->norm,lits);
1309     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1310     /* Test for convergence, xnorm = || X || */
1311     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
1312     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1313     if (snes->reason) break;
1314   }
1315   ierr = DMDestroyVI(snes->dm);CHKERRQ(ierr);
1316   if (i == maxits) {
1317     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
1318     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1319   }
1320   PetscFunctionReturn(0);
1321 }
1322 
1323 #undef __FUNCT__
1324 #define __FUNCT__ "SNESVISetRedundancyCheck"
1325 PetscErrorCode SNESVISetRedundancyCheck(SNES snes,PetscErrorCode (*func)(SNES,IS,IS*,void*),void *ctx)
1326 {
1327   SNES_VI         *vi = (SNES_VI*)snes->data;
1328 
1329   PetscFunctionBegin;
1330   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1331   vi->checkredundancy = func;
1332   vi->ctxP            = ctx;
1333   PetscFunctionReturn(0);
1334 }
1335 
1336 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1337 #include <engine.h>
1338 #include <mex.h>
1339 typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
1340 
1341 #undef __FUNCT__
1342 #define __FUNCT__ "SNESVIRedundancyCheck_Matlab"
1343 PetscErrorCode SNESVIRedundancyCheck_Matlab(SNES snes,IS is_act,IS* is_redact,void* ctx)
1344 {
1345   PetscErrorCode      ierr;
1346   SNESMatlabContext   *sctx = (SNESMatlabContext*)ctx;
1347   int                 nlhs = 1, nrhs = 5;
1348   mxArray             *plhs[1], *prhs[5];
1349   long long int       l1 = 0, l2 = 0, ls = 0;
1350   PetscInt            *indices=PETSC_NULL;
1351 
1352   PetscFunctionBegin;
1353   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1354   PetscValidHeaderSpecific(is_act,IS_CLASSID,2);
1355   PetscValidPointer(is_redact,3);
1356   PetscCheckSameComm(snes,1,is_act,2);
1357 
1358   /* Create IS for reduced active set of size 0, its size and indices will
1359    bet set by the Matlab function */
1360   ierr = ISCreateGeneral(((PetscObject)snes)->comm,0,indices,PETSC_OWN_POINTER,is_redact);CHKERRQ(ierr);
1361   /* call Matlab function in ctx */
1362   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
1363   ierr = PetscMemcpy(&l1,&is_act,sizeof(is_act));CHKERRQ(ierr);
1364   ierr = PetscMemcpy(&l2,is_redact,sizeof(is_act));CHKERRQ(ierr);
1365   prhs[0] = mxCreateDoubleScalar((double)ls);
1366   prhs[1] = mxCreateDoubleScalar((double)l1);
1367   prhs[2] = mxCreateDoubleScalar((double)l2);
1368   prhs[3] = mxCreateString(sctx->funcname);
1369   prhs[4] = sctx->ctx;
1370   ierr    = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESVIRedundancyCheckInternal");CHKERRQ(ierr);
1371   ierr    = mxGetScalar(plhs[0]);CHKERRQ(ierr);
1372   mxDestroyArray(prhs[0]);
1373   mxDestroyArray(prhs[1]);
1374   mxDestroyArray(prhs[2]);
1375   mxDestroyArray(prhs[3]);
1376   mxDestroyArray(plhs[0]);
1377   PetscFunctionReturn(0);
1378 }
1379 
1380 #undef __FUNCT__
1381 #define __FUNCT__ "SNESVISetRedundancyCheckMatlab"
1382 PetscErrorCode SNESVISetRedundancyCheckMatlab(SNES snes,const char* func,mxArray* ctx)
1383 {
1384   PetscErrorCode      ierr;
1385   SNESMatlabContext   *sctx;
1386 
1387   PetscFunctionBegin;
1388   /* currently sctx is memory bleed */
1389   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
1390   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
1391   sctx->ctx = mxDuplicateArray(ctx);
1392   ierr = SNESVISetRedundancyCheck(snes,SNESVIRedundancyCheck_Matlab,sctx);CHKERRQ(ierr);
1393   PetscFunctionReturn(0);
1394 }
1395 
1396 #endif
1397 
1398 
1399 /* Variational Inequality solver using augmented space method. It does the opposite of the
1400    reduced space method i.e. it identifies the active set variables and instead of discarding
1401    them it augments the original system by introducing additional equality
1402    constraint equations for active set variables. The user can optionally provide an IS for
1403    a subset of 'redundant' active set variables which will treated as free variables.
1404    Specific implementation for Allen-Cahn problem
1405 */
1406 #undef __FUNCT__
1407 #define __FUNCT__ "SNESSolveVI_RSAUG"
1408 PetscErrorCode SNESSolveVI_RSAUG(SNES snes)
1409 {
1410   SNES_VI          *vi = (SNES_VI*)snes->data;
1411   PetscErrorCode    ierr;
1412   PetscInt          maxits,i,lits;
1413   PetscBool         lssucceed;
1414   MatStructure      flg = DIFFERENT_NONZERO_PATTERN;
1415   PetscReal         fnorm,gnorm,xnorm=0,ynorm;
1416   Vec                Y,X,F,G,W;
1417   KSPConvergedReason kspreason;
1418 
1419   PetscFunctionBegin;
1420   snes->numFailures            = 0;
1421   snes->numLinearSolveFailures = 0;
1422   snes->reason                 = SNES_CONVERGED_ITERATING;
1423 
1424   maxits	= snes->max_its;	/* maximum number of iterations */
1425   X		= snes->vec_sol;	/* solution vector */
1426   F		= snes->vec_func;	/* residual vector */
1427   Y		= snes->work[0];	/* work vectors */
1428   G		= snes->work[1];
1429   W		= snes->work[2];
1430 
1431   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1432   snes->iter = 0;
1433   snes->norm = 0.0;
1434   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1435 
1436   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
1437   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1438   if (snes->domainerror) {
1439     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1440     PetscFunctionReturn(0);
1441   }
1442   ierr = SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);CHKERRQ(ierr);
1443   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
1444   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
1445   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1446 
1447   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1448   snes->norm = fnorm;
1449   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1450   SNESLogConvHistory(snes,fnorm,0);
1451   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1452 
1453   /* set parameter for default relative tolerance convergence test */
1454   snes->ttol = fnorm*snes->rtol;
1455   /* test convergence */
1456   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1457   if (snes->reason) PetscFunctionReturn(0);
1458 
1459   for (i=0; i<maxits; i++) {
1460     IS             IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
1461     IS             IS_redact; /* redundant active set */
1462     Mat            J_aug,Jpre_aug;
1463     Vec            F_aug,Y_aug;
1464     PetscInt       nis_redact,nis_act;
1465     const PetscInt *idx_redact,*idx_act;
1466     PetscInt       k;
1467     PetscInt       *idx_actkept=PETSC_NULL,nkept=0; /* list of kept active set */
1468     PetscScalar    *f,*f2;
1469     PetscBool      isequal;
1470     PetscInt       i1=0,j1=0;
1471 
1472     /* Call general purpose update function */
1473     if (snes->ops->update) {
1474       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1475     }
1476     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
1477 
1478     /* Create active and inactive index sets */
1479     ierr = SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&IS_inact);CHKERRQ(ierr);
1480 
1481     /* Get local active set size */
1482     ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr);
1483     if (nis_act) {
1484       ierr = ISGetIndices(IS_act,&idx_act);CHKERRQ(ierr);
1485       IS_redact  = PETSC_NULL;
1486       if(vi->checkredundancy) {
1487 	(*vi->checkredundancy)(snes,IS_act,&IS_redact,vi->ctxP);
1488       }
1489 
1490       if(!IS_redact) {
1491 	/* User called checkredundancy function but didn't create IS_redact because
1492            there were no redundant active set variables */
1493 	/* Copy over all active set indices to the list */
1494 	ierr = PetscMalloc(nis_act*sizeof(PetscInt),&idx_actkept);CHKERRQ(ierr);
1495 	for(k=0;k < nis_act;k++) idx_actkept[k] = idx_act[k];
1496 	nkept = nis_act;
1497       } else {
1498 	ierr = ISGetLocalSize(IS_redact,&nis_redact);CHKERRQ(ierr);
1499 	ierr = PetscMalloc((nis_act-nis_redact)*sizeof(PetscInt),&idx_actkept);CHKERRQ(ierr);
1500 
1501 	/* Create reduced active set list */
1502 	ierr = ISGetIndices(IS_act,&idx_act);CHKERRQ(ierr);
1503 	ierr = ISGetIndices(IS_redact,&idx_redact);CHKERRQ(ierr);
1504 	j1 = 0;nkept = 0;
1505 	for(k=0;k<nis_act;k++) {
1506 	  if(j1 < nis_redact && idx_act[k] == idx_redact[j1]) j1++;
1507 	  else idx_actkept[nkept++] = idx_act[k];
1508 	}
1509 	ierr = ISRestoreIndices(IS_act,&idx_act);CHKERRQ(ierr);
1510 	ierr = ISRestoreIndices(IS_redact,&idx_redact);CHKERRQ(ierr);
1511 
1512         ierr = ISDestroy(&IS_redact);CHKERRQ(ierr);
1513       }
1514 
1515       /* Create augmented F and Y */
1516       ierr = VecCreate(((PetscObject)snes)->comm,&F_aug);CHKERRQ(ierr);
1517       ierr = VecSetSizes(F_aug,F->map->n+nkept,PETSC_DECIDE);CHKERRQ(ierr);
1518       ierr = VecSetFromOptions(F_aug);CHKERRQ(ierr);
1519       ierr = VecDuplicate(F_aug,&Y_aug);CHKERRQ(ierr);
1520 
1521       /* Copy over F to F_aug in the first n locations */
1522       ierr = VecGetArray(F,&f);CHKERRQ(ierr);
1523       ierr = VecGetArray(F_aug,&f2);CHKERRQ(ierr);
1524       ierr = PetscMemcpy(f2,f,F->map->n*sizeof(PetscScalar));CHKERRQ(ierr);
1525       ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
1526       ierr = VecRestoreArray(F_aug,&f2);CHKERRQ(ierr);
1527 
1528       /* Create the augmented jacobian matrix */
1529       ierr = MatCreate(((PetscObject)snes)->comm,&J_aug);CHKERRQ(ierr);
1530       ierr = MatSetSizes(J_aug,X->map->n+nkept,X->map->n+nkept,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
1531       ierr = MatSetFromOptions(J_aug);CHKERRQ(ierr);
1532 
1533 
1534       { /* local vars */
1535       /* Preallocate augmented matrix and set values in it...Doing only sequential case first*/
1536       PetscInt          ncols;
1537       const PetscInt    *cols;
1538       const PetscScalar *vals;
1539       PetscScalar        value[2];
1540       PetscInt           row,col[2];
1541       PetscInt           *d_nnz;
1542       value[0] = 1.0; value[1] = 0.0;
1543       ierr = PetscMalloc((X->map->n+nkept)*sizeof(PetscInt),&d_nnz);CHKERRQ(ierr);
1544       ierr = PetscMemzero(d_nnz,(X->map->n+nkept)*sizeof(PetscInt));CHKERRQ(ierr);
1545       for(row=0;row<snes->jacobian->rmap->n;row++) {
1546         ierr = MatGetRow(snes->jacobian,row,&ncols,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1547         d_nnz[row] += ncols;
1548         ierr = MatRestoreRow(snes->jacobian,row,&ncols,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1549       }
1550 
1551       for(k=0;k<nkept;k++) {
1552         d_nnz[idx_actkept[k]] += 1;
1553         d_nnz[snes->jacobian->rmap->n+k] += 2;
1554       }
1555       ierr = MatSeqAIJSetPreallocation(J_aug,PETSC_NULL,d_nnz);CHKERRQ(ierr);
1556 
1557       ierr = PetscFree(d_nnz);CHKERRQ(ierr);
1558 
1559       /* Set the original jacobian matrix in J_aug */
1560       for(row=0;row<snes->jacobian->rmap->n;row++) {
1561 	ierr = MatGetRow(snes->jacobian,row,&ncols,&cols,&vals);CHKERRQ(ierr);
1562 	ierr = MatSetValues(J_aug,1,&row,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
1563 	ierr = MatRestoreRow(snes->jacobian,row,&ncols,&cols,&vals);CHKERRQ(ierr);
1564       }
1565       /* Add the augmented part */
1566       for(k=0;k<nkept;k++) {
1567 	row = snes->jacobian->rmap->n+k;
1568 	col[0] = idx_actkept[k]; col[1] = snes->jacobian->rmap->n+k;
1569 	ierr = MatSetValues(J_aug,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr);
1570 	ierr = MatSetValues(J_aug,1,&col[0],1,&row,&value[0],INSERT_VALUES);CHKERRQ(ierr);
1571       }
1572       ierr = MatAssemblyBegin(J_aug,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1573       ierr = MatAssemblyEnd(J_aug,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1574       /* Only considering prejac = jac for now */
1575       Jpre_aug = J_aug;
1576       } /* local vars*/
1577       ierr = ISRestoreIndices(IS_act,&idx_act);CHKERRQ(ierr);
1578       ierr = ISDestroy(&IS_act);CHKERRQ(ierr);
1579     } else {
1580       F_aug = F; J_aug = snes->jacobian; Y_aug = Y; Jpre_aug = snes->jacobian_pre;
1581     }
1582 
1583     ierr = ISEqual(vi->IS_inact_prev,IS_inact,&isequal);CHKERRQ(ierr);
1584     if (!isequal) {
1585       ierr = SNESVIResetPCandKSP(snes,J_aug,Jpre_aug);CHKERRQ(ierr);
1586     }
1587     ierr = KSPSetOperators(snes->ksp,J_aug,Jpre_aug,flg);CHKERRQ(ierr);
1588     ierr = KSPSetUp(snes->ksp);CHKERRQ(ierr);
1589     /*  {
1590       PC        pc;
1591       PetscBool flg;
1592       ierr = KSPGetPC(snes->ksp,&pc);CHKERRQ(ierr);
1593       ierr = PetscTypeCompare((PetscObject)pc,PCFIELDSPLIT,&flg);CHKERRQ(ierr);
1594       if (flg) {
1595         KSP      *subksps;
1596         ierr = PCFieldSplitGetSubKSP(pc,PETSC_NULL,&subksps);CHKERRQ(ierr);
1597         ierr = KSPGetPC(subksps[0],&pc);CHKERRQ(ierr);
1598         ierr = PetscFree(subksps);CHKERRQ(ierr);
1599         ierr = PetscTypeCompare((PetscObject)pc,PCBJACOBI,&flg);CHKERRQ(ierr);
1600         if (flg) {
1601           PetscInt       n,N = 101*101,j,cnts[3] = {0,0,0};
1602           const PetscInt *ii;
1603 
1604           ierr = ISGetSize(IS_inact,&n);CHKERRQ(ierr);
1605           ierr = ISGetIndices(IS_inact,&ii);CHKERRQ(ierr);
1606           for (j=0; j<n; j++) {
1607             if (ii[j] < N) cnts[0]++;
1608             else if (ii[j] < 2*N) cnts[1]++;
1609             else if (ii[j] < 3*N) cnts[2]++;
1610           }
1611           ierr = ISRestoreIndices(IS_inact,&ii);CHKERRQ(ierr);
1612 
1613           ierr = PCBJacobiSetTotalBlocks(pc,3,cnts);CHKERRQ(ierr);
1614         }
1615       }
1616     }
1617     */
1618     ierr = SNES_KSPSolve(snes,snes->ksp,F_aug,Y_aug);CHKERRQ(ierr);
1619     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
1620     if (kspreason < 0) {
1621       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
1622         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
1623         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
1624         break;
1625       }
1626     }
1627 
1628     if(nis_act) {
1629       PetscScalar *y1,*y2;
1630       ierr = VecGetArray(Y,&y1);CHKERRQ(ierr);
1631       ierr = VecGetArray(Y_aug,&y2);CHKERRQ(ierr);
1632       /* Copy over inactive Y_aug to Y */
1633       j1 = 0;
1634       for(i1=Y->map->rstart;i1 < Y->map->rend;i1++) {
1635 	if(j1 < nkept && idx_actkept[j1] == i1) j1++;
1636 	else y1[i1-Y->map->rstart] = y2[i1-Y->map->rstart];
1637       }
1638       ierr = VecRestoreArray(Y,&y1);CHKERRQ(ierr);
1639       ierr = VecRestoreArray(Y_aug,&y2);CHKERRQ(ierr);
1640 
1641       ierr = VecDestroy(&F_aug);CHKERRQ(ierr);
1642       ierr = VecDestroy(&Y_aug);CHKERRQ(ierr);
1643       ierr = MatDestroy(&J_aug);CHKERRQ(ierr);
1644       ierr = PetscFree(idx_actkept);CHKERRQ(ierr);
1645     }
1646 
1647     if (!isequal) {
1648       ierr = ISDestroy(&vi->IS_inact_prev);CHKERRQ(ierr);
1649       ierr = ISDuplicate(IS_inact,&vi->IS_inact_prev);CHKERRQ(ierr);
1650     }
1651     ierr = ISDestroy(&IS_inact);CHKERRQ(ierr);
1652 
1653     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
1654     snes->linear_its += lits;
1655     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
1656     /*
1657     if (vi->precheckstep) {
1658       PetscBool changed_y = PETSC_FALSE;
1659       ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr);
1660     }
1661 
1662     if (PetscLogPrintInfo){
1663       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
1664     }
1665     */
1666     /* Compute a (scaled) negative update in the line search routine:
1667          Y <- X - lambda*Y
1668        and evaluate G = function(Y) (depends on the line search).
1669     */
1670     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
1671     ynorm = 1; gnorm = fnorm;
1672     ierr = (*vi->LineSearch)(snes,vi->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1673     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
1674     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
1675     if (snes->domainerror) {
1676       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1677       PetscFunctionReturn(0);
1678     }
1679     if (!lssucceed) {
1680       if (++snes->numFailures >= snes->maxFailures) {
1681 	PetscBool ismin;
1682         snes->reason = SNES_DIVERGED_LINE_SEARCH;
1683         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
1684         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
1685         break;
1686       }
1687     }
1688     /* Update function and solution vectors */
1689     fnorm = gnorm;
1690     ierr = VecCopy(G,F);CHKERRQ(ierr);
1691     ierr = VecCopy(W,X);CHKERRQ(ierr);
1692     /* Monitor convergence */
1693     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1694     snes->iter = i+1;
1695     snes->norm = fnorm;
1696     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1697     SNESLogConvHistory(snes,snes->norm,lits);
1698     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1699     /* Test for convergence, xnorm = || X || */
1700     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
1701     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1702     if (snes->reason) break;
1703   }
1704   if (i == maxits) {
1705     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
1706     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1707   }
1708   PetscFunctionReturn(0);
1709 }
1710 
1711 /* -------------------------------------------------------------------------- */
1712 /*
1713    SNESSetUp_VI - Sets up the internal data structures for the later use
1714    of the SNESVI nonlinear solver.
1715 
1716    Input Parameter:
1717 .  snes - the SNES context
1718 .  x - the solution vector
1719 
1720    Application Interface Routine: SNESSetUp()
1721 
1722    Notes:
1723    For basic use of the SNES solvers, the user need not explicitly call
1724    SNESSetUp(), since these actions will automatically occur during
1725    the call to SNESSolve().
1726  */
1727 #undef __FUNCT__
1728 #define __FUNCT__ "SNESSetUp_VI"
1729 PetscErrorCode SNESSetUp_VI(SNES snes)
1730 {
1731   PetscErrorCode ierr;
1732   SNES_VI        *vi = (SNES_VI*) snes->data;
1733   PetscInt       i_start[3],i_end[3];
1734 
1735   PetscFunctionBegin;
1736 
1737   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
1738 
1739   if (vi->computevariablebounds) {
1740     if (!vi->xl) {ierr = VecDuplicate(snes->vec_sol,&vi->xl);CHKERRQ(ierr);}
1741     if (!vi->xu) {ierr = VecDuplicate(snes->vec_sol,&vi->xu);CHKERRQ(ierr);}
1742     ierr = (*vi->computevariablebounds)(snes,vi->xl,vi->xu);CHKERRQ(ierr);
1743   } else if (!vi->xl && !vi->xu) {
1744     /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
1745     ierr = VecDuplicate(snes->vec_sol, &vi->xl);CHKERRQ(ierr);
1746     ierr = VecSet(vi->xl,SNES_VI_NINF);CHKERRQ(ierr);
1747     ierr = VecDuplicate(snes->vec_sol, &vi->xu);CHKERRQ(ierr);
1748     ierr = VecSet(vi->xu,SNES_VI_INF);CHKERRQ(ierr);
1749   } else {
1750     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
1751     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
1752     ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr);
1753     ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr);
1754     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]))
1755       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.");
1756   }
1757   if (snes->ops->solve != SNESSolveVI_SS) {
1758     /* Set up previous active index set for the first snes solve
1759        vi->IS_inact_prev = 0,1,2,....N */
1760     PetscInt *indices;
1761     PetscInt  i,n,rstart,rend;
1762 
1763     ierr = VecGetOwnershipRange(snes->vec_sol,&rstart,&rend);CHKERRQ(ierr);
1764     ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
1765     ierr = PetscMalloc(n*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1766     for(i=0;i < n; i++) indices[i] = rstart + i;
1767     ierr = ISCreateGeneral(((PetscObject)snes)->comm,n,indices,PETSC_OWN_POINTER,&vi->IS_inact_prev);
1768   }
1769 
1770   if (snes->ops->solve == SNESSolveVI_SS) {
1771     ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr);
1772     ierr = VecDuplicate(snes->vec_sol, &vi->phi);CHKERRQ(ierr);
1773     ierr = VecDuplicate(snes->vec_sol, &vi->Da);CHKERRQ(ierr);
1774     ierr = VecDuplicate(snes->vec_sol, &vi->Db);CHKERRQ(ierr);
1775     ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr);
1776     ierr = VecDuplicate(snes->vec_sol, &vi->t);CHKERRQ(ierr);
1777   }
1778   PetscFunctionReturn(0);
1779 }
1780 /* -------------------------------------------------------------------------- */
1781 #undef __FUNCT__
1782 #define __FUNCT__ "SNESReset_VI"
1783 PetscErrorCode SNESReset_VI(SNES snes)
1784 {
1785   SNES_VI        *vi = (SNES_VI*) snes->data;
1786   PetscErrorCode ierr;
1787 
1788   PetscFunctionBegin;
1789   ierr = VecDestroy(&vi->xl);CHKERRQ(ierr);
1790   ierr = VecDestroy(&vi->xu);CHKERRQ(ierr);
1791   if (snes->ops->solve != SNESSolveVI_SS) {
1792     ierr = ISDestroy(&vi->IS_inact_prev);CHKERRQ(ierr);
1793   }
1794   PetscFunctionReturn(0);
1795 }
1796 
1797 /*
1798    SNESDestroy_VI - Destroys the private SNES_VI context that was created
1799    with SNESCreate_VI().
1800 
1801    Input Parameter:
1802 .  snes - the SNES context
1803 
1804    Application Interface Routine: SNESDestroy()
1805  */
1806 #undef __FUNCT__
1807 #define __FUNCT__ "SNESDestroy_VI"
1808 PetscErrorCode SNESDestroy_VI(SNES snes)
1809 {
1810   SNES_VI        *vi = (SNES_VI*) snes->data;
1811   PetscErrorCode ierr;
1812 
1813   PetscFunctionBegin;
1814 
1815   if (snes->ops->solve == SNESSolveVI_SS) {
1816     /* clear vectors */
1817     ierr = VecDestroy(&vi->dpsi);CHKERRQ(ierr);
1818     ierr = VecDestroy(&vi->phi);CHKERRQ(ierr);
1819     ierr = VecDestroy(&vi->Da);CHKERRQ(ierr);
1820     ierr = VecDestroy(&vi->Db);CHKERRQ(ierr);
1821     ierr = VecDestroy(&vi->z);CHKERRQ(ierr);
1822     ierr = VecDestroy(&vi->t);CHKERRQ(ierr);
1823   }
1824 
1825   ierr = PetscViewerDestroy(&vi->lsmonitor);CHKERRQ(ierr);
1826   ierr = PetscFree(snes->data);CHKERRQ(ierr);
1827 
1828   /* clear composed functions */
1829   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
1830   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
1831   PetscFunctionReturn(0);
1832 }
1833 
1834 /* -------------------------------------------------------------------------- */
1835 #undef __FUNCT__
1836 #define __FUNCT__ "SNESLineSearchNo_VI"
1837 
1838 /*
1839   This routine does not actually do a line search but it takes a full newton
1840   step while ensuring that the new iterates remain within the constraints.
1841 
1842 */
1843 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)
1844 {
1845   PetscErrorCode ierr;
1846   SNES_VI        *vi = (SNES_VI*)snes->data;
1847   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1848 
1849   PetscFunctionBegin;
1850   *flag = PETSC_TRUE;
1851   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1852   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
1853   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1854   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1855   if (vi->postcheckstep) {
1856    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1857   }
1858   if (changed_y) {
1859     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1860     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1861   }
1862   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1863   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1864   if (!snes->domainerror) {
1865     if (snes->ops->solve != SNESSolveVI_SS) {
1866        ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1867     } else {
1868       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
1869     }
1870     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1871   }
1872   if (vi->lsmonitor) {
1873     ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1874     ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1875     ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1876   }
1877   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1878   PetscFunctionReturn(0);
1879 }
1880 
1881 /* -------------------------------------------------------------------------- */
1882 #undef __FUNCT__
1883 #define __FUNCT__ "SNESLineSearchNoNorms_VI"
1884 
1885 /*
1886   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
1887 */
1888 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)
1889 {
1890   PetscErrorCode ierr;
1891   SNES_VI        *vi = (SNES_VI*)snes->data;
1892   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1893 
1894   PetscFunctionBegin;
1895   *flag = PETSC_TRUE;
1896   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1897   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
1898   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1899   if (vi->postcheckstep) {
1900    ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
1901   }
1902   if (changed_y) {
1903     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
1904     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1905   }
1906 
1907   /* don't evaluate function the last time through */
1908   if (snes->iter < snes->max_its-1) {
1909     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1910   }
1911   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1912   PetscFunctionReturn(0);
1913 }
1914 
1915 /* -------------------------------------------------------------------------- */
1916 #undef __FUNCT__
1917 #define __FUNCT__ "SNESLineSearchCubic_VI"
1918 /*
1919   This routine implements a cubic line search while doing a projection on the variable bounds
1920 */
1921 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)
1922 {
1923   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
1924   PetscReal      minlambda,lambda,lambdatemp;
1925 #if defined(PETSC_USE_COMPLEX)
1926   PetscScalar    cinitslope;
1927 #endif
1928   PetscErrorCode ierr;
1929   PetscInt       count;
1930   SNES_VI        *vi = (SNES_VI*)snes->data;
1931   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
1932   MPI_Comm       comm;
1933 
1934   PetscFunctionBegin;
1935   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1936   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1937   *flag   = PETSC_TRUE;
1938 
1939   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
1940   if (*ynorm == 0.0) {
1941     if (vi->lsmonitor) {
1942       ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1943       ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
1944       ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1945     }
1946     *gnorm = fnorm;
1947     ierr   = VecCopy(x,w);CHKERRQ(ierr);
1948     ierr   = VecCopy(f,g);CHKERRQ(ierr);
1949     *flag  = PETSC_FALSE;
1950     goto theend1;
1951   }
1952   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
1953     if (vi->lsmonitor) {
1954       ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1955       ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"    Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr);
1956       ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1957     }
1958     ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
1959     *ynorm = vi->maxstep;
1960   }
1961   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
1962   minlambda = vi->minlambda/rellength;
1963   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
1964 #if defined(PETSC_USE_COMPLEX)
1965   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
1966   initslope = PetscRealPart(cinitslope);
1967 #else
1968   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
1969 #endif
1970   if (initslope > 0.0)  initslope = -initslope;
1971   if (initslope == 0.0) initslope = -1.0;
1972 
1973   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
1974   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
1975   if (snes->nfuncs >= snes->max_funcs) {
1976     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
1977     *flag = PETSC_FALSE;
1978     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1979     goto theend1;
1980   }
1981   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
1982   if (snes->ops->solve != SNESSolveVI_SS) {
1983     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
1984   } else {
1985     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
1986   }
1987   if (snes->domainerror) {
1988     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
1989     PetscFunctionReturn(0);
1990   }
1991   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1992   ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1993   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
1994     if (vi->lsmonitor) {
1995       ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1996       ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
1997       ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1998     }
1999     goto theend1;
2000   }
2001 
2002   /* Fit points with quadratic */
2003   lambda     = 1.0;
2004   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
2005   lambdaprev = lambda;
2006   gnormprev  = *gnorm;
2007   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
2008   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
2009   else                         lambda = lambdatemp;
2010 
2011   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
2012   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
2013   if (snes->nfuncs >= snes->max_funcs) {
2014     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
2015     *flag = PETSC_FALSE;
2016     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
2017     goto theend1;
2018   }
2019   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
2020   if (snes->ops->solve != SNESSolveVI_SS) {
2021     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
2022   } else {
2023     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
2024   }
2025   if (snes->domainerror) {
2026     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
2027     PetscFunctionReturn(0);
2028   }
2029   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
2030   if (vi->lsmonitor) {
2031     ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2032     ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"    Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr);
2033     ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2034   }
2035   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
2036     if (vi->lsmonitor) {
2037       ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2038       ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
2039       ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2040     }
2041     goto theend1;
2042   }
2043 
2044   /* Fit points with cubic */
2045   count = 1;
2046   while (PETSC_TRUE) {
2047     if (lambda <= minlambda) {
2048       if (vi->lsmonitor) {
2049         ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2050  	ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
2051 	ierr = PetscViewerASCIIPrintf(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);
2052         ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2053       }
2054       *flag = PETSC_FALSE;
2055       break;
2056     }
2057     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
2058     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
2059     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
2060     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
2061     d  = b*b - 3*a*initslope;
2062     if (d < 0.0) d = 0.0;
2063     if (a == 0.0) {
2064       lambdatemp = -initslope/(2.0*b);
2065     } else {
2066       lambdatemp = (-b + sqrt(d))/(3.0*a);
2067     }
2068     lambdaprev = lambda;
2069     gnormprev  = *gnorm;
2070     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
2071     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
2072     else                         lambda     = lambdatemp;
2073     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
2074     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
2075     if (snes->nfuncs >= snes->max_funcs) {
2076       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
2077       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);
2078       *flag = PETSC_FALSE;
2079       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
2080       break;
2081     }
2082     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
2083     if (snes->ops->solve != SNESSolveVI_SS) {
2084       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
2085     } else {
2086       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
2087     }
2088     if (snes->domainerror) {
2089       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
2090       PetscFunctionReturn(0);
2091     }
2092     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
2093     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */
2094       if (vi->lsmonitor) {
2095 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
2096       }
2097       break;
2098     } else {
2099       if (vi->lsmonitor) {
2100         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
2101       }
2102     }
2103     count++;
2104   }
2105   theend1:
2106   /* Optional user-defined check for line search step validity */
2107   if (vi->postcheckstep && *flag) {
2108     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
2109     if (changed_y) {
2110       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
2111       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
2112     }
2113     if (changed_y || changed_w) { /* recompute the function if the step has changed */
2114       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
2115       if (snes->ops->solve != SNESSolveVI_SS) {
2116         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
2117       } else {
2118         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
2119       }
2120       if (snes->domainerror) {
2121         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
2122         PetscFunctionReturn(0);
2123       }
2124       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
2125       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
2126       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
2127     }
2128   }
2129   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
2130   PetscFunctionReturn(0);
2131 }
2132 
2133 #undef __FUNCT__
2134 #define __FUNCT__ "SNESLineSearchQuadratic_VI"
2135 /*
2136   This routine does a quadratic line search while keeping the iterates within the variable bounds
2137 */
2138 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)
2139 {
2140   /*
2141      Note that for line search purposes we work with with the related
2142      minimization problem:
2143         min  z(x):  R^n -> R,
2144      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
2145    */
2146   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
2147 #if defined(PETSC_USE_COMPLEX)
2148   PetscScalar    cinitslope;
2149 #endif
2150   PetscErrorCode ierr;
2151   PetscInt       count;
2152   SNES_VI        *vi = (SNES_VI*)snes->data;
2153   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
2154 
2155   PetscFunctionBegin;
2156   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
2157   *flag   = PETSC_TRUE;
2158 
2159   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
2160   if (*ynorm == 0.0) {
2161     if (vi->lsmonitor) {
2162       ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2163       ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
2164       ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2165     }
2166     *gnorm = fnorm;
2167     ierr   = VecCopy(x,w);CHKERRQ(ierr);
2168     ierr   = VecCopy(f,g);CHKERRQ(ierr);
2169     *flag  = PETSC_FALSE;
2170     goto theend2;
2171   }
2172   if (*ynorm > vi->maxstep) {	/* Step too big, so scale back */
2173     ierr   = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr);
2174     *ynorm = vi->maxstep;
2175   }
2176   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
2177   minlambda = vi->minlambda/rellength;
2178   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
2179 #if defined(PETSC_USE_COMPLEX)
2180   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
2181   initslope = PetscRealPart(cinitslope);
2182 #else
2183   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
2184 #endif
2185   if (initslope > 0.0)  initslope = -initslope;
2186   if (initslope == 0.0) initslope = -1.0;
2187   ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr);
2188 
2189   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
2190   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
2191   if (snes->nfuncs >= snes->max_funcs) {
2192     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
2193     *flag = PETSC_FALSE;
2194     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
2195     goto theend2;
2196   }
2197   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
2198   if (snes->ops->solve != SNESSolveVI_SS) {
2199     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
2200   } else {
2201     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
2202   }
2203   if (snes->domainerror) {
2204     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
2205     PetscFunctionReturn(0);
2206   }
2207   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
2208   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */
2209     if (vi->lsmonitor) {
2210       ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2211       ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
2212       ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2213     }
2214     goto theend2;
2215   }
2216 
2217   /* Fit points with quadratic */
2218   lambda = 1.0;
2219   count = 1;
2220   while (PETSC_TRUE) {
2221     if (lambda <= minlambda) { /* bad luck; use full step */
2222       if (vi->lsmonitor) {
2223         ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2224         ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
2225         ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"Line search: fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
2226         ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2227       }
2228       ierr = VecCopy(x,w);CHKERRQ(ierr);
2229       *flag = PETSC_FALSE;
2230       break;
2231     }
2232     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
2233     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
2234     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
2235     else                         lambda     = lambdatemp;
2236 
2237     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
2238     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
2239     if (snes->nfuncs >= snes->max_funcs) {
2240       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
2241       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);
2242       *flag = PETSC_FALSE;
2243       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
2244       break;
2245     }
2246     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
2247     if (snes->domainerror) {
2248       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
2249       PetscFunctionReturn(0);
2250     }
2251     if (snes->ops->solve != SNESSolveVI_SS) {
2252       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
2253     } else {
2254       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
2255     }
2256     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
2257     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */
2258       if (vi->lsmonitor) {
2259         ierr = PetscViewerASCIIAddTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2260         ierr = PetscViewerASCIIPrintf(vi->lsmonitor,"    Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
2261         ierr = PetscViewerASCIISubtractTab(vi->lsmonitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
2262       }
2263       break;
2264     }
2265     count++;
2266   }
2267   theend2:
2268   /* Optional user-defined check for line search step validity */
2269   if (vi->postcheckstep) {
2270     ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
2271     if (changed_y) {
2272       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
2273       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
2274     }
2275     if (changed_y || changed_w) { /* recompute the function if the step has changed */
2276       ierr = SNESComputeFunction(snes,w,g);
2277       if (snes->domainerror) {
2278         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
2279         PetscFunctionReturn(0);
2280       }
2281       if (snes->ops->solve != SNESSolveVI_SS) {
2282         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
2283       } else {
2284         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
2285       }
2286 
2287       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
2288       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
2289       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
2290     }
2291   }
2292   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
2293   PetscFunctionReturn(0);
2294 }
2295 
2296 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*/
2297 /* -------------------------------------------------------------------------- */
2298 EXTERN_C_BEGIN
2299 #undef __FUNCT__
2300 #define __FUNCT__ "SNESLineSearchSet_VI"
2301 PetscErrorCode  SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx)
2302 {
2303   PetscFunctionBegin;
2304   ((SNES_VI *)(snes->data))->LineSearch = func;
2305   ((SNES_VI *)(snes->data))->lsP        = lsctx;
2306   PetscFunctionReturn(0);
2307 }
2308 EXTERN_C_END
2309 
2310 /* -------------------------------------------------------------------------- */
2311 EXTERN_C_BEGIN
2312 #undef __FUNCT__
2313 #define __FUNCT__ "SNESLineSearchSetMonitor_VI"
2314 PetscErrorCode  SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg)
2315 {
2316   SNES_VI        *vi = (SNES_VI*)snes->data;
2317   PetscErrorCode ierr;
2318 
2319   PetscFunctionBegin;
2320   if (flg && !vi->lsmonitor) {
2321     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&vi->lsmonitor);CHKERRQ(ierr);
2322   } else if (!flg && vi->lsmonitor) {
2323     ierr = PetscViewerDestroy(&vi->lsmonitor);CHKERRQ(ierr);
2324   }
2325   PetscFunctionReturn(0);
2326 }
2327 EXTERN_C_END
2328 
2329 /*
2330    SNESView_VI - Prints info from the SNESVI data structure.
2331 
2332    Input Parameters:
2333 .  SNES - the SNES context
2334 .  viewer - visualization context
2335 
2336    Application Interface Routine: SNESView()
2337 */
2338 #undef __FUNCT__
2339 #define __FUNCT__ "SNESView_VI"
2340 static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer)
2341 {
2342   SNES_VI        *vi = (SNES_VI *)snes->data;
2343   const char     *cstr,*tstr;
2344   PetscErrorCode ierr;
2345   PetscBool     iascii;
2346 
2347   PetscFunctionBegin;
2348   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2349   if (iascii) {
2350     if (vi->LineSearch == SNESLineSearchNo_VI)             cstr = "SNESLineSearchNo";
2351     else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic";
2352     else if (vi->LineSearch == SNESLineSearchCubic_VI)     cstr = "SNESLineSearchCubic";
2353     else                                                   cstr = "unknown";
2354     if (snes->ops->solve == SNESSolveVI_SS)         tstr = "Semismooth";
2355     else if (snes->ops->solve == SNESSolveVI_RS)    tstr = "Reduced Space";
2356     else if (snes->ops->solve == SNESSolveVI_RSAUG) tstr = "Reduced space with augmented variables";
2357     else                                         tstr = "unknown";
2358     ierr = PetscViewerASCIIPrintf(viewer,"  VI algorithm: %s\n",tstr);CHKERRQ(ierr);
2359     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
2360     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr);
2361   } else {
2362     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name);
2363   }
2364   PetscFunctionReturn(0);
2365 }
2366 
2367 #undef __FUNCT__
2368 #define __FUNCT__ "SNESVISetVariableBounds"
2369 /*@
2370    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
2371 
2372    Input Parameters:
2373 .  snes - the SNES context.
2374 .  xl   - lower bound.
2375 .  xu   - upper bound.
2376 
2377    Notes:
2378    If this routine is not called then the lower and upper bounds are set to
2379    SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp().
2380 
2381 @*/
2382 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
2383 {
2384   SNES_VI        *vi;
2385   PetscErrorCode ierr;
2386 
2387   PetscFunctionBegin;
2388   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2389   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
2390   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
2391   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
2392   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);
2393   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);
2394   ierr = SNESSetType(snes,SNESVI);CHKERRQ(ierr);
2395   vi = (SNES_VI*)snes->data;
2396   ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
2397   ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
2398   ierr = VecDestroy(&vi->xl);CHKERRQ(ierr);
2399   ierr = VecDestroy(&vi->xu);CHKERRQ(ierr);
2400   vi->xl = xl;
2401   vi->xu = xu;
2402   PetscFunctionReturn(0);
2403 }
2404 
2405 /* -------------------------------------------------------------------------- */
2406 /*
2407    SNESSetFromOptions_VI - Sets various parameters for the SNESVI method.
2408 
2409    Input Parameter:
2410 .  snes - the SNES context
2411 
2412    Application Interface Routine: SNESSetFromOptions()
2413 */
2414 #undef __FUNCT__
2415 #define __FUNCT__ "SNESSetFromOptions_VI"
2416 static PetscErrorCode SNESSetFromOptions_VI(SNES snes)
2417 {
2418   SNES_VI        *vi = (SNES_VI *)snes->data;
2419   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
2420   const char     *vies[] = {"ss","rs","rsaug"};
2421   PetscErrorCode ierr;
2422   PetscInt       indx;
2423   PetscBool      flg,set,flg2;
2424 
2425   PetscFunctionBegin;
2426   ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr);
2427   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
2428   if (flg) {
2429     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
2430   }
2431   ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr);
2432   ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr);
2433   ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr);
2434   ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr);
2435   ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
2436   if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
2437   ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,3,"ss",&indx,&flg2);CHKERRQ(ierr);
2438   if (flg2) {
2439     switch (indx) {
2440     case 0:
2441       snes->ops->solve = SNESSolveVI_SS;
2442       break;
2443     case 1:
2444       snes->ops->solve = SNESSolveVI_RS;
2445       break;
2446     case 2:
2447       snes->ops->solve = SNESSolveVI_RSAUG;
2448     }
2449   }
2450   ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"basic",&indx,&flg);CHKERRQ(ierr);
2451   if (flg) {
2452     switch (indx) {
2453     case 0:
2454       ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
2455       break;
2456     case 1:
2457       ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
2458       break;
2459     case 2:
2460       ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
2461       break;
2462     case 3:
2463       ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
2464       break;
2465     }
2466   }
2467   ierr = PetscOptionsTail();CHKERRQ(ierr);
2468   PetscFunctionReturn(0);
2469 }
2470 /* -------------------------------------------------------------------------- */
2471 /*MC
2472       SNESVI - Various solvers for variational inequalities based on Newton's method
2473 
2474    Options Database:
2475 +   -snes_vi_type <ss,rs,rsaug> a semi-smooth solver, a reduced space active set method and a reduced space active set method that does not eliminate the active constraints from the Jacobian instead augments the Jacobian with
2476                                 additional variables that enforce the constraints
2477 -   -snes_vi_monitor - prints the number of active constraints at each iteration.
2478 
2479 
2480    Level: beginner
2481 
2482 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
2483            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
2484            SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
2485 
2486 M*/
2487 EXTERN_C_BEGIN
2488 #undef __FUNCT__
2489 #define __FUNCT__ "SNESCreate_VI"
2490 PetscErrorCode  SNESCreate_VI(SNES snes)
2491 {
2492   PetscErrorCode ierr;
2493   SNES_VI        *vi;
2494 
2495   PetscFunctionBegin;
2496   snes->ops->reset           = SNESReset_VI;
2497   snes->ops->setup           = SNESSetUp_VI;
2498   snes->ops->solve           = SNESSolveVI_RS;
2499   snes->ops->destroy         = SNESDestroy_VI;
2500   snes->ops->setfromoptions  = SNESSetFromOptions_VI;
2501   snes->ops->view            = SNESView_VI;
2502   snes->ops->converged       = SNESDefaultConverged_VI;
2503 
2504   ierr                  = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr);
2505   snes->data            = (void*)vi;
2506   vi->alpha             = 1.e-4;
2507   vi->maxstep           = 1.e8;
2508   vi->minlambda         = 1.e-12;
2509   vi->LineSearch        = SNESLineSearchNo_VI;
2510   vi->lsP               = PETSC_NULL;
2511   vi->postcheckstep     = PETSC_NULL;
2512   vi->postcheck         = PETSC_NULL;
2513   vi->precheckstep      = PETSC_NULL;
2514   vi->precheck          = PETSC_NULL;
2515   vi->const_tol         =  2.2204460492503131e-16;
2516   vi->checkredundancy   = PETSC_NULL;
2517 
2518   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr);
2519   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr);
2520 
2521   PetscFunctionReturn(0);
2522 }
2523 EXTERN_C_END
2524