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