xref: /petsc/src/snes/impls/vi/rs/virs.c (revision c2fc9fa91a764ea1bc445c5db6fc55248852cdef)
1*c2fc9fa9SBarry Smith 
2*c2fc9fa9SBarry Smith #include <../src/snes/impls/vi/rs/virsimpl.h> /*I "petscsnes.h" I*/
3*c2fc9fa9SBarry Smith #include <../include/private/kspimpl.h>
4*c2fc9fa9SBarry Smith #include <../include/private/matimpl.h>
5*c2fc9fa9SBarry Smith #include <../include/private/dmimpl.h>
6*c2fc9fa9SBarry Smith 
7*c2fc9fa9SBarry Smith #undef __FUNCT__
8*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIGetInactiveSet"
9*c2fc9fa9SBarry Smith /*
10*c2fc9fa9SBarry Smith    SNESVIGetInactiveSet - Gets the global indices for the inactive set variables (these correspond to the degrees of freedom the linear
11*c2fc9fa9SBarry Smith      system is solved on)
12*c2fc9fa9SBarry Smith 
13*c2fc9fa9SBarry Smith    Input parameter
14*c2fc9fa9SBarry Smith .  snes - the SNES context
15*c2fc9fa9SBarry Smith 
16*c2fc9fa9SBarry Smith    Output parameter
17*c2fc9fa9SBarry Smith .  ISact - active set index set
18*c2fc9fa9SBarry Smith 
19*c2fc9fa9SBarry Smith  */
20*c2fc9fa9SBarry Smith PetscErrorCode SNESVIGetInactiveSet(SNES snes,IS* inact)
21*c2fc9fa9SBarry Smith {
22*c2fc9fa9SBarry Smith   SNES_VIRS        *vi = (SNES_VIRS*)snes->data;
23*c2fc9fa9SBarry Smith   PetscFunctionBegin;
24*c2fc9fa9SBarry Smith   *inact = vi->IS_inact_prev;
25*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
26*c2fc9fa9SBarry Smith }
27*c2fc9fa9SBarry Smith 
28*c2fc9fa9SBarry Smith /*
29*c2fc9fa9SBarry Smith     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
30*c2fc9fa9SBarry Smith   defined by the reduced space method.
31*c2fc9fa9SBarry Smith 
32*c2fc9fa9SBarry Smith     Simple calls the regular DM interpolation and restricts it to operation on the variables not associated with active constraints.
33*c2fc9fa9SBarry Smith 
34*c2fc9fa9SBarry Smith <*/
35*c2fc9fa9SBarry Smith typedef struct {
36*c2fc9fa9SBarry Smith   PetscInt       n;                                        /* size of vectors in the reduced DM space */
37*c2fc9fa9SBarry Smith   IS             inactive;
38*c2fc9fa9SBarry Smith   PetscErrorCode (*getinterpolation)(DM,DM,Mat*,Vec*);    /* DM's original routines */
39*c2fc9fa9SBarry Smith   PetscErrorCode (*coarsen)(DM, MPI_Comm, DM*);
40*c2fc9fa9SBarry Smith   PetscErrorCode (*createglobalvector)(DM,Vec*);
41*c2fc9fa9SBarry Smith   DM             dm;                                      /* when destroying this object we need to reset the above function into the base DM */
42*c2fc9fa9SBarry Smith } DM_SNESVI;
43*c2fc9fa9SBarry Smith 
44*c2fc9fa9SBarry Smith #undef __FUNCT__
45*c2fc9fa9SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_SNESVI"
46*c2fc9fa9SBarry Smith /*
47*c2fc9fa9SBarry Smith      DMCreateGlobalVector_SNESVI - Creates global vector of the size of the reduced space
48*c2fc9fa9SBarry Smith 
49*c2fc9fa9SBarry Smith */
50*c2fc9fa9SBarry Smith PetscErrorCode  DMCreateGlobalVector_SNESVI(DM dm,Vec *vec)
51*c2fc9fa9SBarry Smith {
52*c2fc9fa9SBarry Smith   PetscErrorCode          ierr;
53*c2fc9fa9SBarry Smith   PetscContainer          isnes;
54*c2fc9fa9SBarry Smith   DM_SNESVI               *dmsnesvi;
55*c2fc9fa9SBarry Smith 
56*c2fc9fa9SBarry Smith   PetscFunctionBegin;
57*c2fc9fa9SBarry Smith   ierr = PetscObjectQuery((PetscObject)dm,"VI",(PetscObject *)&isnes);CHKERRQ(ierr);
58*c2fc9fa9SBarry Smith   if (!isnes) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_PLIB,"Composed SNES is missing");
59*c2fc9fa9SBarry Smith   ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi);CHKERRQ(ierr);
60*c2fc9fa9SBarry Smith   ierr = VecCreateMPI(((PetscObject)dm)->comm,dmsnesvi->n,PETSC_DETERMINE,vec);CHKERRQ(ierr);
61*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
62*c2fc9fa9SBarry Smith }
63*c2fc9fa9SBarry Smith 
64*c2fc9fa9SBarry Smith #undef __FUNCT__
65*c2fc9fa9SBarry Smith #define __FUNCT__ "DMGetInterpolation_SNESVI"
66*c2fc9fa9SBarry Smith /*
67*c2fc9fa9SBarry Smith      DMGetInterpolation_SNESVI - Modifieds the interpolation obtained from the DM by removing all rows and columns associated with active constraints.
68*c2fc9fa9SBarry Smith 
69*c2fc9fa9SBarry Smith */
70*c2fc9fa9SBarry Smith PetscErrorCode  DMGetInterpolation_SNESVI(DM dm1,DM dm2,Mat *mat,Vec *vec)
71*c2fc9fa9SBarry Smith {
72*c2fc9fa9SBarry Smith   PetscErrorCode          ierr;
73*c2fc9fa9SBarry Smith   PetscContainer          isnes;
74*c2fc9fa9SBarry Smith   DM_SNESVI               *dmsnesvi1,*dmsnesvi2;
75*c2fc9fa9SBarry Smith   Mat                     interp;
76*c2fc9fa9SBarry Smith 
77*c2fc9fa9SBarry Smith   PetscFunctionBegin;
78*c2fc9fa9SBarry Smith   ierr = PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject *)&isnes);CHKERRQ(ierr);
79*c2fc9fa9SBarry Smith   if (!isnes) SETERRQ(((PetscObject)dm1)->comm,PETSC_ERR_PLIB,"Composed VI data structure is missing");
80*c2fc9fa9SBarry Smith   ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);CHKERRQ(ierr);
81*c2fc9fa9SBarry Smith   ierr = PetscObjectQuery((PetscObject)dm2,"VI",(PetscObject *)&isnes);CHKERRQ(ierr);
82*c2fc9fa9SBarry Smith   if (!isnes) SETERRQ(((PetscObject)dm2)->comm,PETSC_ERR_PLIB,"Composed VI data structure is missing");
83*c2fc9fa9SBarry Smith   ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi2);CHKERRQ(ierr);
84*c2fc9fa9SBarry Smith 
85*c2fc9fa9SBarry Smith   ierr = (*dmsnesvi1->getinterpolation)(dm1,dm2,&interp,PETSC_NULL);CHKERRQ(ierr);
86*c2fc9fa9SBarry Smith   ierr = MatGetSubMatrix(interp,dmsnesvi2->inactive,dmsnesvi1->inactive,MAT_INITIAL_MATRIX,mat);CHKERRQ(ierr);
87*c2fc9fa9SBarry Smith   ierr = MatDestroy(&interp);CHKERRQ(ierr);
88*c2fc9fa9SBarry Smith   *vec = 0;
89*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
90*c2fc9fa9SBarry Smith }
91*c2fc9fa9SBarry Smith 
92*c2fc9fa9SBarry Smith extern PetscErrorCode  DMSetVI(DM,IS);
93*c2fc9fa9SBarry Smith 
94*c2fc9fa9SBarry Smith #undef __FUNCT__
95*c2fc9fa9SBarry Smith #define __FUNCT__ "DMCoarsen_SNESVI"
96*c2fc9fa9SBarry Smith /*
97*c2fc9fa9SBarry Smith      DMCoarsen_SNESVI - Computes the regular coarsened DM then computes additional information about its inactive set
98*c2fc9fa9SBarry Smith 
99*c2fc9fa9SBarry Smith */
100*c2fc9fa9SBarry Smith PetscErrorCode  DMCoarsen_SNESVI(DM dm1,MPI_Comm comm,DM *dm2)
101*c2fc9fa9SBarry Smith {
102*c2fc9fa9SBarry Smith   PetscErrorCode          ierr;
103*c2fc9fa9SBarry Smith   PetscContainer          isnes;
104*c2fc9fa9SBarry Smith   DM_SNESVI               *dmsnesvi1;
105*c2fc9fa9SBarry Smith   Vec                     finemarked,coarsemarked;
106*c2fc9fa9SBarry Smith   IS                      inactive;
107*c2fc9fa9SBarry Smith   VecScatter              inject;
108*c2fc9fa9SBarry Smith   const PetscInt          *index;
109*c2fc9fa9SBarry Smith   PetscInt                n,k,cnt = 0,rstart,*coarseindex;
110*c2fc9fa9SBarry Smith   PetscScalar             *marked;
111*c2fc9fa9SBarry Smith 
112*c2fc9fa9SBarry Smith   PetscFunctionBegin;
113*c2fc9fa9SBarry Smith   ierr = PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject *)&isnes);CHKERRQ(ierr);
114*c2fc9fa9SBarry Smith   if (!isnes) SETERRQ(((PetscObject)dm1)->comm,PETSC_ERR_PLIB,"Composed VI data structure is missing");
115*c2fc9fa9SBarry Smith   ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);CHKERRQ(ierr);
116*c2fc9fa9SBarry Smith 
117*c2fc9fa9SBarry Smith   /* get the original coarsen */
118*c2fc9fa9SBarry Smith   ierr = (*dmsnesvi1->coarsen)(dm1,comm,dm2);CHKERRQ(ierr);
119*c2fc9fa9SBarry Smith 
120*c2fc9fa9SBarry Smith   /* not sure why this extra reference is needed, but without the dm2 disappears too early */
121*c2fc9fa9SBarry Smith   ierr = PetscObjectReference((PetscObject)*dm2);CHKERRQ(ierr);
122*c2fc9fa9SBarry Smith 
123*c2fc9fa9SBarry Smith   /* need to set back global vectors in order to use the original injection */
124*c2fc9fa9SBarry Smith   ierr = DMClearGlobalVectors(dm1);CHKERRQ(ierr);
125*c2fc9fa9SBarry Smith   dm1->ops->createglobalvector = dmsnesvi1->createglobalvector;
126*c2fc9fa9SBarry Smith   ierr = DMCreateGlobalVector(dm1,&finemarked);CHKERRQ(ierr);
127*c2fc9fa9SBarry Smith   ierr = DMCreateGlobalVector(*dm2,&coarsemarked);CHKERRQ(ierr);
128*c2fc9fa9SBarry Smith 
129*c2fc9fa9SBarry Smith   /*
130*c2fc9fa9SBarry Smith      fill finemarked with locations of inactive points
131*c2fc9fa9SBarry Smith   */
132*c2fc9fa9SBarry Smith   ierr = ISGetIndices(dmsnesvi1->inactive,&index);CHKERRQ(ierr);
133*c2fc9fa9SBarry Smith   ierr = ISGetLocalSize(dmsnesvi1->inactive,&n);CHKERRQ(ierr);
134*c2fc9fa9SBarry Smith   ierr = VecSet(finemarked,0.0);CHKERRQ(ierr);
135*c2fc9fa9SBarry Smith   for (k=0;k<n;k++){
136*c2fc9fa9SBarry Smith       ierr = VecSetValue(finemarked,index[k],1.0,INSERT_VALUES);CHKERRQ(ierr);
137*c2fc9fa9SBarry Smith   }
138*c2fc9fa9SBarry Smith   ierr = VecAssemblyBegin(finemarked);CHKERRQ(ierr);
139*c2fc9fa9SBarry Smith   ierr = VecAssemblyEnd(finemarked);CHKERRQ(ierr);
140*c2fc9fa9SBarry Smith 
141*c2fc9fa9SBarry Smith   ierr = DMGetInjection(*dm2,dm1,&inject);CHKERRQ(ierr);
142*c2fc9fa9SBarry Smith   ierr = VecScatterBegin(inject,finemarked,coarsemarked,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
143*c2fc9fa9SBarry Smith   ierr = VecScatterEnd(inject,finemarked,coarsemarked,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
144*c2fc9fa9SBarry Smith   ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
145*c2fc9fa9SBarry Smith 
146*c2fc9fa9SBarry Smith   /*
147*c2fc9fa9SBarry Smith      create index set list of coarse inactive points from coarsemarked
148*c2fc9fa9SBarry Smith   */
149*c2fc9fa9SBarry Smith   ierr = VecGetLocalSize(coarsemarked,&n);CHKERRQ(ierr);
150*c2fc9fa9SBarry Smith   ierr = VecGetOwnershipRange(coarsemarked,&rstart,PETSC_NULL);CHKERRQ(ierr);
151*c2fc9fa9SBarry Smith   ierr = VecGetArray(coarsemarked,&marked);CHKERRQ(ierr);
152*c2fc9fa9SBarry Smith   for (k=0; k<n; k++) {
153*c2fc9fa9SBarry Smith     if (marked[k] != 0.0) cnt++;
154*c2fc9fa9SBarry Smith   }
155*c2fc9fa9SBarry Smith   ierr = PetscMalloc(cnt*sizeof(PetscInt),&coarseindex);CHKERRQ(ierr);
156*c2fc9fa9SBarry Smith   cnt  = 0;
157*c2fc9fa9SBarry Smith   for (k=0; k<n; k++) {
158*c2fc9fa9SBarry Smith     if (marked[k] != 0.0) coarseindex[cnt++] = k + rstart;
159*c2fc9fa9SBarry Smith   }
160*c2fc9fa9SBarry Smith   ierr = VecRestoreArray(coarsemarked,&marked);CHKERRQ(ierr);
161*c2fc9fa9SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_WORLD,cnt,coarseindex,PETSC_OWN_POINTER,&inactive);CHKERRQ(ierr);
162*c2fc9fa9SBarry Smith 
163*c2fc9fa9SBarry Smith   ierr = DMClearGlobalVectors(dm1);CHKERRQ(ierr);
164*c2fc9fa9SBarry Smith   dm1->ops->createglobalvector = DMCreateGlobalVector_SNESVI;
165*c2fc9fa9SBarry Smith   ierr = DMSetVI(*dm2,inactive);CHKERRQ(ierr);
166*c2fc9fa9SBarry Smith 
167*c2fc9fa9SBarry Smith   ierr = VecDestroy(&finemarked);CHKERRQ(ierr);
168*c2fc9fa9SBarry Smith   ierr = VecDestroy(&coarsemarked);CHKERRQ(ierr);
169*c2fc9fa9SBarry Smith   ierr = ISDestroy(&inactive);CHKERRQ(ierr);
170*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
171*c2fc9fa9SBarry Smith }
172*c2fc9fa9SBarry Smith 
173*c2fc9fa9SBarry Smith #undef __FUNCT__
174*c2fc9fa9SBarry Smith #define __FUNCT__ "DMDestroy_SNESVI"
175*c2fc9fa9SBarry Smith PetscErrorCode DMDestroy_SNESVI(DM_SNESVI *dmsnesvi)
176*c2fc9fa9SBarry Smith {
177*c2fc9fa9SBarry Smith   PetscErrorCode ierr;
178*c2fc9fa9SBarry Smith 
179*c2fc9fa9SBarry Smith   PetscFunctionBegin;
180*c2fc9fa9SBarry Smith   /* reset the base methods in the DM object that were changed when the DM_SNESVI was reset */
181*c2fc9fa9SBarry Smith   dmsnesvi->dm->ops->getinterpolation   = dmsnesvi->getinterpolation;
182*c2fc9fa9SBarry Smith   dmsnesvi->dm->ops->coarsen            = dmsnesvi->coarsen;
183*c2fc9fa9SBarry Smith   dmsnesvi->dm->ops->createglobalvector = dmsnesvi->createglobalvector;
184*c2fc9fa9SBarry Smith   /* need to clear out this vectors because some of them may not have a reference to the DM
185*c2fc9fa9SBarry Smith     but they are counted as having references to the DM in DMDestroy() */
186*c2fc9fa9SBarry Smith   ierr = DMClearGlobalVectors(dmsnesvi->dm);CHKERRQ(ierr);
187*c2fc9fa9SBarry Smith 
188*c2fc9fa9SBarry Smith   ierr = ISDestroy(&dmsnesvi->inactive);CHKERRQ(ierr);
189*c2fc9fa9SBarry Smith   ierr = PetscFree(dmsnesvi);CHKERRQ(ierr);
190*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
191*c2fc9fa9SBarry Smith }
192*c2fc9fa9SBarry Smith 
193*c2fc9fa9SBarry Smith #undef __FUNCT__
194*c2fc9fa9SBarry Smith #define __FUNCT__ "DMSetVI"
195*c2fc9fa9SBarry Smith /*
196*c2fc9fa9SBarry Smith      DMSetVI - Marks a DM as associated with a VI problem. This causes the interpolation/restriction operators to
197*c2fc9fa9SBarry Smith                be restricted to only those variables NOT associated with active constraints.
198*c2fc9fa9SBarry Smith 
199*c2fc9fa9SBarry Smith */
200*c2fc9fa9SBarry Smith PetscErrorCode  DMSetVI(DM dm,IS inactive)
201*c2fc9fa9SBarry Smith {
202*c2fc9fa9SBarry Smith   PetscErrorCode          ierr;
203*c2fc9fa9SBarry Smith   PetscContainer          isnes;
204*c2fc9fa9SBarry Smith   DM_SNESVI               *dmsnesvi;
205*c2fc9fa9SBarry Smith 
206*c2fc9fa9SBarry Smith   PetscFunctionBegin;
207*c2fc9fa9SBarry Smith   if (!dm) PetscFunctionReturn(0);
208*c2fc9fa9SBarry Smith 
209*c2fc9fa9SBarry Smith   ierr = PetscObjectReference((PetscObject)inactive);CHKERRQ(ierr);
210*c2fc9fa9SBarry Smith 
211*c2fc9fa9SBarry Smith   ierr = PetscObjectQuery((PetscObject)dm,"VI",(PetscObject *)&isnes);CHKERRQ(ierr);
212*c2fc9fa9SBarry Smith   if (!isnes) {
213*c2fc9fa9SBarry Smith     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&isnes);CHKERRQ(ierr);
214*c2fc9fa9SBarry Smith     ierr = PetscContainerSetUserDestroy(isnes,(PetscErrorCode (*)(void*))DMDestroy_SNESVI);CHKERRQ(ierr);
215*c2fc9fa9SBarry Smith     ierr = PetscNew(DM_SNESVI,&dmsnesvi);CHKERRQ(ierr);
216*c2fc9fa9SBarry Smith     ierr = PetscContainerSetPointer(isnes,(void*)dmsnesvi);CHKERRQ(ierr);
217*c2fc9fa9SBarry Smith     ierr = PetscObjectCompose((PetscObject)dm,"VI",(PetscObject)isnes);CHKERRQ(ierr);
218*c2fc9fa9SBarry Smith     ierr = PetscContainerDestroy(&isnes);CHKERRQ(ierr);
219*c2fc9fa9SBarry Smith     dmsnesvi->getinterpolation   = dm->ops->getinterpolation;
220*c2fc9fa9SBarry Smith     dm->ops->getinterpolation    = DMGetInterpolation_SNESVI;
221*c2fc9fa9SBarry Smith     dmsnesvi->coarsen            = dm->ops->coarsen;
222*c2fc9fa9SBarry Smith     dm->ops->coarsen             = DMCoarsen_SNESVI;
223*c2fc9fa9SBarry Smith     dmsnesvi->createglobalvector = dm->ops->createglobalvector;
224*c2fc9fa9SBarry Smith     dm->ops->createglobalvector  = DMCreateGlobalVector_SNESVI;
225*c2fc9fa9SBarry Smith   } else {
226*c2fc9fa9SBarry Smith     ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi);CHKERRQ(ierr);
227*c2fc9fa9SBarry Smith     ierr = ISDestroy(&dmsnesvi->inactive);CHKERRQ(ierr);
228*c2fc9fa9SBarry Smith   }
229*c2fc9fa9SBarry Smith   ierr = DMClearGlobalVectors(dm);CHKERRQ(ierr);
230*c2fc9fa9SBarry Smith   ierr = ISGetLocalSize(inactive,&dmsnesvi->n);CHKERRQ(ierr);
231*c2fc9fa9SBarry Smith   dmsnesvi->inactive = inactive;
232*c2fc9fa9SBarry Smith   dmsnesvi->dm       = dm;
233*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
234*c2fc9fa9SBarry Smith }
235*c2fc9fa9SBarry Smith 
236*c2fc9fa9SBarry Smith #undef __FUNCT__
237*c2fc9fa9SBarry Smith #define __FUNCT__ "DMDestroyVI"
238*c2fc9fa9SBarry Smith /*
239*c2fc9fa9SBarry Smith      DMDestroyVI - Frees the DM_SNESVI object contained in the DM
240*c2fc9fa9SBarry Smith          - also resets the function pointers in the DM for getinterpolation() etc to use the original DM
241*c2fc9fa9SBarry Smith */
242*c2fc9fa9SBarry Smith PetscErrorCode  DMDestroyVI(DM dm)
243*c2fc9fa9SBarry Smith {
244*c2fc9fa9SBarry Smith   PetscErrorCode          ierr;
245*c2fc9fa9SBarry Smith 
246*c2fc9fa9SBarry Smith   PetscFunctionBegin;
247*c2fc9fa9SBarry Smith   if (!dm) PetscFunctionReturn(0);
248*c2fc9fa9SBarry Smith   ierr = PetscObjectCompose((PetscObject)dm,"VI",(PetscObject)PETSC_NULL);CHKERRQ(ierr);
249*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
250*c2fc9fa9SBarry Smith }
251*c2fc9fa9SBarry Smith 
252*c2fc9fa9SBarry Smith /* --------------------------------------------------------------------------------------------------------*/
253*c2fc9fa9SBarry Smith 
254*c2fc9fa9SBarry Smith 
255*c2fc9fa9SBarry Smith 
256*c2fc9fa9SBarry Smith 
257*c2fc9fa9SBarry Smith #undef __FUNCT__
258*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESCreateIndexSets_VIRS"
259*c2fc9fa9SBarry Smith PetscErrorCode SNESCreateIndexSets_VIRS(SNES snes,Vec X,Vec F,IS* ISact,IS* ISinact)
260*c2fc9fa9SBarry Smith {
261*c2fc9fa9SBarry Smith   PetscErrorCode     ierr;
262*c2fc9fa9SBarry Smith 
263*c2fc9fa9SBarry Smith   PetscFunctionBegin;
264*c2fc9fa9SBarry Smith   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
265*c2fc9fa9SBarry Smith   ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr);
266*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
267*c2fc9fa9SBarry Smith }
268*c2fc9fa9SBarry Smith 
269*c2fc9fa9SBarry Smith /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */
270*c2fc9fa9SBarry Smith #undef __FUNCT__
271*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESCreateSubVectors_VIRS"
272*c2fc9fa9SBarry Smith PetscErrorCode SNESCreateSubVectors_VIRS(SNES snes,PetscInt n,Vec* newv)
273*c2fc9fa9SBarry Smith {
274*c2fc9fa9SBarry Smith   PetscErrorCode ierr;
275*c2fc9fa9SBarry Smith   Vec            v;
276*c2fc9fa9SBarry Smith 
277*c2fc9fa9SBarry Smith   PetscFunctionBegin;
278*c2fc9fa9SBarry Smith   ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr);
279*c2fc9fa9SBarry Smith   ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr);
280*c2fc9fa9SBarry Smith   ierr = VecSetFromOptions(v);CHKERRQ(ierr);
281*c2fc9fa9SBarry Smith   *newv = v;
282*c2fc9fa9SBarry Smith 
283*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
284*c2fc9fa9SBarry Smith }
285*c2fc9fa9SBarry Smith 
286*c2fc9fa9SBarry Smith /* Resets the snes PC and KSP when the active set sizes change */
287*c2fc9fa9SBarry Smith #undef __FUNCT__
288*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIResetPCandKSP"
289*c2fc9fa9SBarry Smith PetscErrorCode SNESVIResetPCandKSP(SNES snes,Mat Amat,Mat Pmat)
290*c2fc9fa9SBarry Smith {
291*c2fc9fa9SBarry Smith   PetscErrorCode         ierr;
292*c2fc9fa9SBarry Smith   KSP                    snesksp;
293*c2fc9fa9SBarry Smith 
294*c2fc9fa9SBarry Smith   PetscFunctionBegin;
295*c2fc9fa9SBarry Smith   ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr);
296*c2fc9fa9SBarry Smith   ierr = KSPReset(snesksp);CHKERRQ(ierr);
297*c2fc9fa9SBarry Smith 
298*c2fc9fa9SBarry Smith   /*
299*c2fc9fa9SBarry Smith   KSP                    kspnew;
300*c2fc9fa9SBarry Smith   PC                     pcnew;
301*c2fc9fa9SBarry Smith   const MatSolverPackage stype;
302*c2fc9fa9SBarry Smith 
303*c2fc9fa9SBarry Smith 
304*c2fc9fa9SBarry Smith   ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr);
305*c2fc9fa9SBarry Smith   kspnew->pc_side = snesksp->pc_side;
306*c2fc9fa9SBarry Smith   kspnew->rtol    = snesksp->rtol;
307*c2fc9fa9SBarry Smith   kspnew->abstol    = snesksp->abstol;
308*c2fc9fa9SBarry Smith   kspnew->max_it  = snesksp->max_it;
309*c2fc9fa9SBarry Smith   ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr);
310*c2fc9fa9SBarry Smith   ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr);
311*c2fc9fa9SBarry Smith   ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr);
312*c2fc9fa9SBarry Smith   ierr = PCSetOperators(kspnew->pc,Amat,Pmat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
313*c2fc9fa9SBarry Smith   ierr = PCFactorGetMatSolverPackage(snesksp->pc,&stype);CHKERRQ(ierr);
314*c2fc9fa9SBarry Smith   ierr = PCFactorSetMatSolverPackage(kspnew->pc,stype);CHKERRQ(ierr);
315*c2fc9fa9SBarry Smith   ierr = KSPDestroy(&snesksp);CHKERRQ(ierr);
316*c2fc9fa9SBarry Smith   snes->ksp = kspnew;
317*c2fc9fa9SBarry Smith   ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr);
318*c2fc9fa9SBarry Smith    ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr);*/
319*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
320*c2fc9fa9SBarry Smith }
321*c2fc9fa9SBarry Smith 
322*c2fc9fa9SBarry Smith 
323*c2fc9fa9SBarry Smith 
324*c2fc9fa9SBarry Smith /* Variational Inequality solver using reduce space method. No semismooth algorithm is
325*c2fc9fa9SBarry Smith    implemented in this algorithm. It basically identifies the active constraints and does
326*c2fc9fa9SBarry Smith    a linear solve on the other variables (those not associated with the active constraints). */
327*c2fc9fa9SBarry Smith #undef __FUNCT__
328*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSolve_VIRS"
329*c2fc9fa9SBarry Smith PetscErrorCode SNESSolve_VIRS(SNES snes)
330*c2fc9fa9SBarry Smith {
331*c2fc9fa9SBarry Smith   SNES_VIRS         *vi = (SNES_VIRS*)snes->data;
332*c2fc9fa9SBarry Smith   PetscErrorCode    ierr;
333*c2fc9fa9SBarry Smith   PetscInt          maxits,i,lits;
334*c2fc9fa9SBarry Smith   PetscBool         lssucceed;
335*c2fc9fa9SBarry Smith   MatStructure      flg = DIFFERENT_NONZERO_PATTERN;
336*c2fc9fa9SBarry Smith   PetscReal         fnorm,gnorm,xnorm=0,ynorm;
337*c2fc9fa9SBarry Smith   Vec                Y,X,F,G,W;
338*c2fc9fa9SBarry Smith   KSPConvergedReason kspreason;
339*c2fc9fa9SBarry Smith 
340*c2fc9fa9SBarry Smith   PetscFunctionBegin;
341*c2fc9fa9SBarry Smith   snes->numFailures            = 0;
342*c2fc9fa9SBarry Smith   snes->numLinearSolveFailures = 0;
343*c2fc9fa9SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
344*c2fc9fa9SBarry Smith 
345*c2fc9fa9SBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
346*c2fc9fa9SBarry Smith   X		= snes->vec_sol;	/* solution vector */
347*c2fc9fa9SBarry Smith   F		= snes->vec_func;	/* residual vector */
348*c2fc9fa9SBarry Smith   Y		= snes->work[0];	/* work vectors */
349*c2fc9fa9SBarry Smith   G		= snes->work[1];
350*c2fc9fa9SBarry Smith   W		= snes->work[2];
351*c2fc9fa9SBarry Smith 
352*c2fc9fa9SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
353*c2fc9fa9SBarry Smith   snes->iter = 0;
354*c2fc9fa9SBarry Smith   snes->norm = 0.0;
355*c2fc9fa9SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
356*c2fc9fa9SBarry Smith 
357*c2fc9fa9SBarry Smith   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
358*c2fc9fa9SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
359*c2fc9fa9SBarry Smith   if (snes->domainerror) {
360*c2fc9fa9SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
361*c2fc9fa9SBarry Smith     PetscFunctionReturn(0);
362*c2fc9fa9SBarry Smith   }
363*c2fc9fa9SBarry Smith   ierr = SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);CHKERRQ(ierr);
364*c2fc9fa9SBarry Smith   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
365*c2fc9fa9SBarry Smith   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
366*c2fc9fa9SBarry Smith   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
367*c2fc9fa9SBarry Smith 
368*c2fc9fa9SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
369*c2fc9fa9SBarry Smith   snes->norm = fnorm;
370*c2fc9fa9SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
371*c2fc9fa9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
372*c2fc9fa9SBarry Smith   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
373*c2fc9fa9SBarry Smith 
374*c2fc9fa9SBarry Smith   /* set parameter for default relative tolerance convergence test */
375*c2fc9fa9SBarry Smith   snes->ttol = fnorm*snes->rtol;
376*c2fc9fa9SBarry Smith   /* test convergence */
377*c2fc9fa9SBarry Smith   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
378*c2fc9fa9SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
379*c2fc9fa9SBarry Smith 
380*c2fc9fa9SBarry Smith 
381*c2fc9fa9SBarry Smith   for (i=0; i<maxits; i++) {
382*c2fc9fa9SBarry Smith 
383*c2fc9fa9SBarry Smith     IS         IS_act,IS_inact; /* _act -> active set _inact -> inactive set */
384*c2fc9fa9SBarry Smith     IS         IS_redact; /* redundant active set */
385*c2fc9fa9SBarry Smith     VecScatter scat_act,scat_inact;
386*c2fc9fa9SBarry Smith     PetscInt   nis_act,nis_inact;
387*c2fc9fa9SBarry Smith     Vec        Y_act,Y_inact,F_inact;
388*c2fc9fa9SBarry Smith     Mat        jac_inact_inact,prejac_inact_inact;
389*c2fc9fa9SBarry Smith     IS         keptrows;
390*c2fc9fa9SBarry Smith     PetscBool  isequal;
391*c2fc9fa9SBarry Smith 
392*c2fc9fa9SBarry Smith     /* Call general purpose update function */
393*c2fc9fa9SBarry Smith     if (snes->ops->update) {
394*c2fc9fa9SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
395*c2fc9fa9SBarry Smith     }
396*c2fc9fa9SBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
397*c2fc9fa9SBarry Smith 
398*c2fc9fa9SBarry Smith 
399*c2fc9fa9SBarry Smith         /* Create active and inactive index sets */
400*c2fc9fa9SBarry Smith 
401*c2fc9fa9SBarry Smith     /*original
402*c2fc9fa9SBarry Smith     ierr = SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&IS_inact);CHKERRQ(ierr);
403*c2fc9fa9SBarry Smith      */
404*c2fc9fa9SBarry Smith     ierr = SNESVIGetActiveSetIS(snes,X,F,&IS_act);CHKERRQ(ierr);
405*c2fc9fa9SBarry Smith 
406*c2fc9fa9SBarry Smith     if (vi->checkredundancy) {
407*c2fc9fa9SBarry Smith       (*vi->checkredundancy)(snes,IS_act,&IS_redact,vi->ctxP);CHKERRQ(ierr);
408*c2fc9fa9SBarry Smith       if (IS_redact){
409*c2fc9fa9SBarry Smith         ierr = ISSort(IS_redact);CHKERRQ(ierr);
410*c2fc9fa9SBarry Smith         ierr = ISComplement(IS_redact,X->map->rstart,X->map->rend,&IS_inact);CHKERRQ(ierr);
411*c2fc9fa9SBarry Smith         ierr = ISDestroy(&IS_redact);CHKERRQ(ierr);
412*c2fc9fa9SBarry Smith       }
413*c2fc9fa9SBarry Smith       else {
414*c2fc9fa9SBarry Smith         ierr = ISComplement(IS_act,X->map->rstart,X->map->rend,&IS_inact);CHKERRQ(ierr);
415*c2fc9fa9SBarry Smith       }
416*c2fc9fa9SBarry Smith     } else {
417*c2fc9fa9SBarry Smith       ierr = ISComplement(IS_act,X->map->rstart,X->map->rend,&IS_inact);CHKERRQ(ierr);
418*c2fc9fa9SBarry Smith     }
419*c2fc9fa9SBarry Smith 
420*c2fc9fa9SBarry Smith 
421*c2fc9fa9SBarry Smith     /* Create inactive set submatrix */
422*c2fc9fa9SBarry Smith     ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
423*c2fc9fa9SBarry Smith 
424*c2fc9fa9SBarry Smith     ierr = MatFindNonzeroRows(jac_inact_inact,&keptrows);CHKERRQ(ierr);
425*c2fc9fa9SBarry Smith     if (0 && keptrows) {
426*c2fc9fa9SBarry Smith       PetscInt       cnt,*nrows,k;
427*c2fc9fa9SBarry Smith       const PetscInt *krows,*inact;
428*c2fc9fa9SBarry Smith       PetscInt       rstart=jac_inact_inact->rmap->rstart;
429*c2fc9fa9SBarry Smith 
430*c2fc9fa9SBarry Smith       ierr = MatDestroy(&jac_inact_inact);CHKERRQ(ierr);
431*c2fc9fa9SBarry Smith       ierr = ISDestroy(&IS_act);CHKERRQ(ierr);
432*c2fc9fa9SBarry Smith 
433*c2fc9fa9SBarry Smith       ierr = ISGetLocalSize(keptrows,&cnt);CHKERRQ(ierr);
434*c2fc9fa9SBarry Smith       ierr = ISGetIndices(keptrows,&krows);CHKERRQ(ierr);
435*c2fc9fa9SBarry Smith       ierr = ISGetIndices(IS_inact,&inact);CHKERRQ(ierr);
436*c2fc9fa9SBarry Smith       ierr = PetscMalloc(cnt*sizeof(PetscInt),&nrows);CHKERRQ(ierr);
437*c2fc9fa9SBarry Smith       for (k=0; k<cnt; k++) {
438*c2fc9fa9SBarry Smith         nrows[k] = inact[krows[k]-rstart];
439*c2fc9fa9SBarry Smith       }
440*c2fc9fa9SBarry Smith       ierr = ISRestoreIndices(keptrows,&krows);CHKERRQ(ierr);
441*c2fc9fa9SBarry Smith       ierr = ISRestoreIndices(IS_inact,&inact);CHKERRQ(ierr);
442*c2fc9fa9SBarry Smith       ierr = ISDestroy(&keptrows);CHKERRQ(ierr);
443*c2fc9fa9SBarry Smith       ierr = ISDestroy(&IS_inact);CHKERRQ(ierr);
444*c2fc9fa9SBarry Smith 
445*c2fc9fa9SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_WORLD,cnt,nrows,PETSC_OWN_POINTER,&IS_inact);CHKERRQ(ierr);
446*c2fc9fa9SBarry Smith       ierr = ISComplement(IS_inact,F->map->rstart,F->map->rend,&IS_act);CHKERRQ(ierr);
447*c2fc9fa9SBarry Smith       ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr);
448*c2fc9fa9SBarry Smith     }
449*c2fc9fa9SBarry Smith     ierr = DMSetVI(snes->dm,IS_inact);CHKERRQ(ierr);
450*c2fc9fa9SBarry Smith     /* remove later */
451*c2fc9fa9SBarry Smith 
452*c2fc9fa9SBarry Smith     /*
453*c2fc9fa9SBarry Smith   ierr = VecView(vi->xu,PETSC_VIEWER_BINARY_(PETSC_COMM_WORLD));CHKERRQ(ierr);
454*c2fc9fa9SBarry Smith   ierr = VecView(vi->xl,PETSC_VIEWER_BINARY_(PETSC_COMM_WORLD));CHKERRQ(ierr);
455*c2fc9fa9SBarry Smith   ierr = VecView(X,PETSC_VIEWER_BINARY_(PETSC_COMM_WORLD));CHKERRQ(ierr);
456*c2fc9fa9SBarry Smith   ierr = VecView(F,PETSC_VIEWER_BINARY_(PETSC_COMM_WORLD));CHKERRQ(ierr);
457*c2fc9fa9SBarry Smith   ierr = ISView(IS_inact,PETSC_VIEWER_BINARY_(PETSC_COMM_WORLD));CHKERRQ(ierr);
458*c2fc9fa9SBarry Smith      */
459*c2fc9fa9SBarry Smith 
460*c2fc9fa9SBarry Smith     /* Get sizes of active and inactive sets */
461*c2fc9fa9SBarry Smith     ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr);
462*c2fc9fa9SBarry Smith     ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr);
463*c2fc9fa9SBarry Smith 
464*c2fc9fa9SBarry Smith     /* Create active and inactive set vectors */
465*c2fc9fa9SBarry Smith     ierr = SNESCreateSubVectors_VIRS(snes,nis_inact,&F_inact);CHKERRQ(ierr);
466*c2fc9fa9SBarry Smith     ierr = SNESCreateSubVectors_VIRS(snes,nis_act,&Y_act);CHKERRQ(ierr);
467*c2fc9fa9SBarry Smith     ierr = SNESCreateSubVectors_VIRS(snes,nis_inact,&Y_inact);CHKERRQ(ierr);
468*c2fc9fa9SBarry Smith 
469*c2fc9fa9SBarry Smith     /* Create scatter contexts */
470*c2fc9fa9SBarry Smith     ierr = VecScatterCreate(Y,IS_act,Y_act,PETSC_NULL,&scat_act);CHKERRQ(ierr);
471*c2fc9fa9SBarry Smith     ierr = VecScatterCreate(Y,IS_inact,Y_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr);
472*c2fc9fa9SBarry Smith 
473*c2fc9fa9SBarry Smith     /* Do a vec scatter to active and inactive set vectors */
474*c2fc9fa9SBarry Smith     ierr = VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
475*c2fc9fa9SBarry Smith     ierr = VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
476*c2fc9fa9SBarry Smith 
477*c2fc9fa9SBarry Smith     ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
478*c2fc9fa9SBarry Smith     ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
479*c2fc9fa9SBarry Smith 
480*c2fc9fa9SBarry Smith     ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
481*c2fc9fa9SBarry Smith     ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
482*c2fc9fa9SBarry Smith 
483*c2fc9fa9SBarry Smith     /* Active set direction = 0 */
484*c2fc9fa9SBarry Smith     ierr = VecSet(Y_act,0);CHKERRQ(ierr);
485*c2fc9fa9SBarry Smith     if (snes->jacobian != snes->jacobian_pre) {
486*c2fc9fa9SBarry Smith       ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr);
487*c2fc9fa9SBarry Smith     } else prejac_inact_inact = jac_inact_inact;
488*c2fc9fa9SBarry Smith 
489*c2fc9fa9SBarry Smith     ierr = ISEqual(vi->IS_inact_prev,IS_inact,&isequal);CHKERRQ(ierr);
490*c2fc9fa9SBarry Smith     if (!isequal) {
491*c2fc9fa9SBarry Smith       ierr = SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);CHKERRQ(ierr);
492*c2fc9fa9SBarry Smith       flg  = DIFFERENT_NONZERO_PATTERN;
493*c2fc9fa9SBarry Smith     }
494*c2fc9fa9SBarry Smith 
495*c2fc9fa9SBarry Smith     /*      ierr = ISView(IS_inact,0);CHKERRQ(ierr); */
496*c2fc9fa9SBarry Smith     /*      ierr = ISView(IS_act,0);CHKERRQ(ierr);*/
497*c2fc9fa9SBarry Smith     /*      ierr = MatView(snes->jacobian_pre,0); */
498*c2fc9fa9SBarry Smith 
499*c2fc9fa9SBarry Smith 
500*c2fc9fa9SBarry Smith 
501*c2fc9fa9SBarry Smith     ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr);
502*c2fc9fa9SBarry Smith     ierr = KSPSetUp(snes->ksp);CHKERRQ(ierr);
503*c2fc9fa9SBarry Smith     {
504*c2fc9fa9SBarry Smith       PC        pc;
505*c2fc9fa9SBarry Smith       PetscBool flg;
506*c2fc9fa9SBarry Smith       ierr = KSPGetPC(snes->ksp,&pc);CHKERRQ(ierr);
507*c2fc9fa9SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc,PCFIELDSPLIT,&flg);CHKERRQ(ierr);
508*c2fc9fa9SBarry Smith       if (flg) {
509*c2fc9fa9SBarry Smith         KSP      *subksps;
510*c2fc9fa9SBarry Smith         ierr = PCFieldSplitGetSubKSP(pc,PETSC_NULL,&subksps);CHKERRQ(ierr);
511*c2fc9fa9SBarry Smith         ierr = KSPGetPC(subksps[0],&pc);CHKERRQ(ierr);
512*c2fc9fa9SBarry Smith         ierr = PetscFree(subksps);CHKERRQ(ierr);
513*c2fc9fa9SBarry Smith         ierr = PetscTypeCompare((PetscObject)pc,PCBJACOBI,&flg);CHKERRQ(ierr);
514*c2fc9fa9SBarry Smith         if (flg) {
515*c2fc9fa9SBarry Smith           PetscInt       n,N = 101*101,j,cnts[3] = {0,0,0};
516*c2fc9fa9SBarry Smith           const PetscInt *ii;
517*c2fc9fa9SBarry Smith 
518*c2fc9fa9SBarry Smith           ierr = ISGetSize(IS_inact,&n);CHKERRQ(ierr);
519*c2fc9fa9SBarry Smith           ierr = ISGetIndices(IS_inact,&ii);CHKERRQ(ierr);
520*c2fc9fa9SBarry Smith           for (j=0; j<n; j++) {
521*c2fc9fa9SBarry Smith             if (ii[j] < N) cnts[0]++;
522*c2fc9fa9SBarry Smith             else if (ii[j] < 2*N) cnts[1]++;
523*c2fc9fa9SBarry Smith             else if (ii[j] < 3*N) cnts[2]++;
524*c2fc9fa9SBarry Smith           }
525*c2fc9fa9SBarry Smith           ierr = ISRestoreIndices(IS_inact,&ii);CHKERRQ(ierr);
526*c2fc9fa9SBarry Smith 
527*c2fc9fa9SBarry Smith           ierr = PCBJacobiSetTotalBlocks(pc,3,cnts);CHKERRQ(ierr);
528*c2fc9fa9SBarry Smith         }
529*c2fc9fa9SBarry Smith       }
530*c2fc9fa9SBarry Smith     }
531*c2fc9fa9SBarry Smith 
532*c2fc9fa9SBarry Smith     ierr = SNES_KSPSolve(snes,snes->ksp,F_inact,Y_inact);CHKERRQ(ierr);
533*c2fc9fa9SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
534*c2fc9fa9SBarry Smith     if (kspreason < 0) {
535*c2fc9fa9SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
536*c2fc9fa9SBarry Smith         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
537*c2fc9fa9SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
538*c2fc9fa9SBarry Smith         break;
539*c2fc9fa9SBarry Smith       }
540*c2fc9fa9SBarry Smith      }
541*c2fc9fa9SBarry Smith 
542*c2fc9fa9SBarry Smith     ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
543*c2fc9fa9SBarry Smith     ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
544*c2fc9fa9SBarry Smith     ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
545*c2fc9fa9SBarry Smith     ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
546*c2fc9fa9SBarry Smith 
547*c2fc9fa9SBarry Smith     ierr = VecDestroy(&F_inact);CHKERRQ(ierr);
548*c2fc9fa9SBarry Smith     ierr = VecDestroy(&Y_act);CHKERRQ(ierr);
549*c2fc9fa9SBarry Smith     ierr = VecDestroy(&Y_inact);CHKERRQ(ierr);
550*c2fc9fa9SBarry Smith     ierr = VecScatterDestroy(&scat_act);CHKERRQ(ierr);
551*c2fc9fa9SBarry Smith     ierr = VecScatterDestroy(&scat_inact);CHKERRQ(ierr);
552*c2fc9fa9SBarry Smith     ierr = ISDestroy(&IS_act);CHKERRQ(ierr);
553*c2fc9fa9SBarry Smith     if (!isequal) {
554*c2fc9fa9SBarry Smith       ierr = ISDestroy(&vi->IS_inact_prev);CHKERRQ(ierr);
555*c2fc9fa9SBarry Smith       ierr = ISDuplicate(IS_inact,&vi->IS_inact_prev);CHKERRQ(ierr);
556*c2fc9fa9SBarry Smith     }
557*c2fc9fa9SBarry Smith     ierr = ISDestroy(&IS_inact);CHKERRQ(ierr);
558*c2fc9fa9SBarry Smith     ierr = MatDestroy(&jac_inact_inact);CHKERRQ(ierr);
559*c2fc9fa9SBarry Smith     if (snes->jacobian != snes->jacobian_pre) {
560*c2fc9fa9SBarry Smith       ierr = MatDestroy(&prejac_inact_inact);CHKERRQ(ierr);
561*c2fc9fa9SBarry Smith     }
562*c2fc9fa9SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
563*c2fc9fa9SBarry Smith     snes->linear_its += lits;
564*c2fc9fa9SBarry Smith     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
565*c2fc9fa9SBarry Smith     /*
566*c2fc9fa9SBarry Smith     if (snes->ops->precheckstep) {
567*c2fc9fa9SBarry Smith       PetscBool changed_y = PETSC_FALSE;
568*c2fc9fa9SBarry Smith       ierr = (*snes->ops->precheckstep)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr);
569*c2fc9fa9SBarry Smith     }
570*c2fc9fa9SBarry Smith 
571*c2fc9fa9SBarry Smith     if (PetscLogPrintInfo){
572*c2fc9fa9SBarry Smith       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
573*c2fc9fa9SBarry Smith     }
574*c2fc9fa9SBarry Smith     */
575*c2fc9fa9SBarry Smith     /* Compute a (scaled) negative update in the line search routine:
576*c2fc9fa9SBarry Smith          Y <- X - lambda*Y
577*c2fc9fa9SBarry Smith        and evaluate G = function(Y) (depends on the line search).
578*c2fc9fa9SBarry Smith     */
579*c2fc9fa9SBarry Smith     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
580*c2fc9fa9SBarry Smith     ynorm = 1; gnorm = fnorm;
581*c2fc9fa9SBarry Smith     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
582*c2fc9fa9SBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
583*c2fc9fa9SBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
584*c2fc9fa9SBarry Smith     if (snes->domainerror) {
585*c2fc9fa9SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
586*c2fc9fa9SBarry Smith       ierr = DMDestroyVI(snes->dm);CHKERRQ(ierr);
587*c2fc9fa9SBarry Smith       PetscFunctionReturn(0);
588*c2fc9fa9SBarry Smith     }
589*c2fc9fa9SBarry Smith     if (!lssucceed) {
590*c2fc9fa9SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
591*c2fc9fa9SBarry Smith 	PetscBool ismin;
592*c2fc9fa9SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
593*c2fc9fa9SBarry Smith         ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
594*c2fc9fa9SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
595*c2fc9fa9SBarry Smith         break;
596*c2fc9fa9SBarry Smith       }
597*c2fc9fa9SBarry Smith     }
598*c2fc9fa9SBarry Smith     /* Update function and solution vectors */
599*c2fc9fa9SBarry Smith     fnorm = gnorm;
600*c2fc9fa9SBarry Smith     ierr = VecCopy(G,F);CHKERRQ(ierr);
601*c2fc9fa9SBarry Smith     ierr = VecCopy(W,X);CHKERRQ(ierr);
602*c2fc9fa9SBarry Smith     /* Monitor convergence */
603*c2fc9fa9SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
604*c2fc9fa9SBarry Smith     snes->iter = i+1;
605*c2fc9fa9SBarry Smith     snes->norm = fnorm;
606*c2fc9fa9SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
607*c2fc9fa9SBarry Smith     SNESLogConvHistory(snes,snes->norm,lits);
608*c2fc9fa9SBarry Smith     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
609*c2fc9fa9SBarry Smith     /* Test for convergence, xnorm = || X || */
610*c2fc9fa9SBarry Smith     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
611*c2fc9fa9SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
612*c2fc9fa9SBarry Smith     if (snes->reason) break;
613*c2fc9fa9SBarry Smith   }
614*c2fc9fa9SBarry Smith   ierr = DMDestroyVI(snes->dm);CHKERRQ(ierr);
615*c2fc9fa9SBarry Smith   if (i == maxits) {
616*c2fc9fa9SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
617*c2fc9fa9SBarry Smith     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
618*c2fc9fa9SBarry Smith   }
619*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
620*c2fc9fa9SBarry Smith }
621*c2fc9fa9SBarry Smith 
622*c2fc9fa9SBarry Smith #undef __FUNCT__
623*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVISetRedundancyCheck"
624*c2fc9fa9SBarry Smith PetscErrorCode SNESVISetRedundancyCheck(SNES snes,PetscErrorCode (*func)(SNES,IS,IS*,void*),void *ctx)
625*c2fc9fa9SBarry Smith {
626*c2fc9fa9SBarry Smith   SNES_VIRS  *vi = (SNES_VIRS*)snes->data;
627*c2fc9fa9SBarry Smith 
628*c2fc9fa9SBarry Smith   PetscFunctionBegin;
629*c2fc9fa9SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
630*c2fc9fa9SBarry Smith   vi->checkredundancy = func;
631*c2fc9fa9SBarry Smith   vi->ctxP            = ctx;
632*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
633*c2fc9fa9SBarry Smith }
634*c2fc9fa9SBarry Smith 
635*c2fc9fa9SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
636*c2fc9fa9SBarry Smith #include <engine.h>
637*c2fc9fa9SBarry Smith #include <mex.h>
638*c2fc9fa9SBarry Smith typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
639*c2fc9fa9SBarry Smith 
640*c2fc9fa9SBarry Smith #undef __FUNCT__
641*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIRedundancyCheck_Matlab"
642*c2fc9fa9SBarry Smith PetscErrorCode SNESVIRedundancyCheck_Matlab(SNES snes,IS is_act,IS* is_redact,void* ctx)
643*c2fc9fa9SBarry Smith {
644*c2fc9fa9SBarry Smith   PetscErrorCode      ierr;
645*c2fc9fa9SBarry Smith   SNESMatlabContext   *sctx = (SNESMatlabContext*)ctx;
646*c2fc9fa9SBarry Smith   int                 nlhs = 1, nrhs = 5;
647*c2fc9fa9SBarry Smith   mxArray             *plhs[1], *prhs[5];
648*c2fc9fa9SBarry Smith   long long int       l1 = 0, l2 = 0, ls = 0;
649*c2fc9fa9SBarry Smith   PetscInt            *indices=PETSC_NULL;
650*c2fc9fa9SBarry Smith 
651*c2fc9fa9SBarry Smith   PetscFunctionBegin;
652*c2fc9fa9SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
653*c2fc9fa9SBarry Smith   PetscValidHeaderSpecific(is_act,IS_CLASSID,2);
654*c2fc9fa9SBarry Smith   PetscValidPointer(is_redact,3);
655*c2fc9fa9SBarry Smith   PetscCheckSameComm(snes,1,is_act,2);
656*c2fc9fa9SBarry Smith 
657*c2fc9fa9SBarry Smith   /* Create IS for reduced active set of size 0, its size and indices will
658*c2fc9fa9SBarry Smith    bet set by the Matlab function */
659*c2fc9fa9SBarry Smith   ierr = ISCreateGeneral(((PetscObject)snes)->comm,0,indices,PETSC_OWN_POINTER,is_redact);CHKERRQ(ierr);
660*c2fc9fa9SBarry Smith   /* call Matlab function in ctx */
661*c2fc9fa9SBarry Smith   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
662*c2fc9fa9SBarry Smith   ierr = PetscMemcpy(&l1,&is_act,sizeof(is_act));CHKERRQ(ierr);
663*c2fc9fa9SBarry Smith   ierr = PetscMemcpy(&l2,is_redact,sizeof(is_act));CHKERRQ(ierr);
664*c2fc9fa9SBarry Smith   prhs[0] = mxCreateDoubleScalar((double)ls);
665*c2fc9fa9SBarry Smith   prhs[1] = mxCreateDoubleScalar((double)l1);
666*c2fc9fa9SBarry Smith   prhs[2] = mxCreateDoubleScalar((double)l2);
667*c2fc9fa9SBarry Smith   prhs[3] = mxCreateString(sctx->funcname);
668*c2fc9fa9SBarry Smith   prhs[4] = sctx->ctx;
669*c2fc9fa9SBarry Smith   ierr    = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESVIRedundancyCheckInternal");CHKERRQ(ierr);
670*c2fc9fa9SBarry Smith   ierr    = mxGetScalar(plhs[0]);CHKERRQ(ierr);
671*c2fc9fa9SBarry Smith   mxDestroyArray(prhs[0]);
672*c2fc9fa9SBarry Smith   mxDestroyArray(prhs[1]);
673*c2fc9fa9SBarry Smith   mxDestroyArray(prhs[2]);
674*c2fc9fa9SBarry Smith   mxDestroyArray(prhs[3]);
675*c2fc9fa9SBarry Smith   mxDestroyArray(plhs[0]);
676*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
677*c2fc9fa9SBarry Smith }
678*c2fc9fa9SBarry Smith 
679*c2fc9fa9SBarry Smith #undef __FUNCT__
680*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVISetRedundancyCheckMatlab"
681*c2fc9fa9SBarry Smith PetscErrorCode SNESVISetRedundancyCheckMatlab(SNES snes,const char* func,mxArray* ctx)
682*c2fc9fa9SBarry Smith {
683*c2fc9fa9SBarry Smith   PetscErrorCode      ierr;
684*c2fc9fa9SBarry Smith   SNESMatlabContext   *sctx;
685*c2fc9fa9SBarry Smith 
686*c2fc9fa9SBarry Smith   PetscFunctionBegin;
687*c2fc9fa9SBarry Smith   /* currently sctx is memory bleed */
688*c2fc9fa9SBarry Smith   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
689*c2fc9fa9SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
690*c2fc9fa9SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
691*c2fc9fa9SBarry Smith   ierr = SNESVISetRedundancyCheck(snes,SNESVIRedundancyCheck_Matlab,sctx);CHKERRQ(ierr);
692*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
693*c2fc9fa9SBarry Smith }
694*c2fc9fa9SBarry Smith 
695*c2fc9fa9SBarry Smith #endif
696*c2fc9fa9SBarry Smith 
697*c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
698*c2fc9fa9SBarry Smith /*
699*c2fc9fa9SBarry Smith    SNESSetUp_VIRS - Sets up the internal data structures for the later use
700*c2fc9fa9SBarry Smith    of the SNESVI nonlinear solver.
701*c2fc9fa9SBarry Smith 
702*c2fc9fa9SBarry Smith    Input Parameter:
703*c2fc9fa9SBarry Smith .  snes - the SNES context
704*c2fc9fa9SBarry Smith .  x - the solution vector
705*c2fc9fa9SBarry Smith 
706*c2fc9fa9SBarry Smith    Application Interface Routine: SNESSetUp()
707*c2fc9fa9SBarry Smith 
708*c2fc9fa9SBarry Smith    Notes:
709*c2fc9fa9SBarry Smith    For basic use of the SNES solvers, the user need not explicitly call
710*c2fc9fa9SBarry Smith    SNESSetUp(), since these actions will automatically occur during
711*c2fc9fa9SBarry Smith    the call to SNESSolve().
712*c2fc9fa9SBarry Smith  */
713*c2fc9fa9SBarry Smith #undef __FUNCT__
714*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetUp_VIRS"
715*c2fc9fa9SBarry Smith PetscErrorCode SNESSetUp_VIRS(SNES snes)
716*c2fc9fa9SBarry Smith {
717*c2fc9fa9SBarry Smith   PetscErrorCode ierr;
718*c2fc9fa9SBarry Smith   SNES_VIRS      *vi = (SNES_VIRS*) snes->data;
719*c2fc9fa9SBarry Smith   PetscInt       *indices;
720*c2fc9fa9SBarry Smith   PetscInt       i,n,rstart,rend;
721*c2fc9fa9SBarry Smith 
722*c2fc9fa9SBarry Smith   PetscFunctionBegin;
723*c2fc9fa9SBarry Smith   ierr = SNESSetUp_VI(snes);CHKERRQ(ierr);
724*c2fc9fa9SBarry Smith 
725*c2fc9fa9SBarry Smith   /* Set up previous active index set for the first snes solve
726*c2fc9fa9SBarry Smith    vi->IS_inact_prev = 0,1,2,....N */
727*c2fc9fa9SBarry Smith 
728*c2fc9fa9SBarry Smith   ierr = VecGetOwnershipRange(snes->vec_sol,&rstart,&rend);CHKERRQ(ierr);
729*c2fc9fa9SBarry Smith   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
730*c2fc9fa9SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&indices);CHKERRQ(ierr);
731*c2fc9fa9SBarry Smith   for(i=0;i < n; i++) indices[i] = rstart + i;
732*c2fc9fa9SBarry Smith   ierr = ISCreateGeneral(((PetscObject)snes)->comm,n,indices,PETSC_OWN_POINTER,&vi->IS_inact_prev);
733*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
734*c2fc9fa9SBarry Smith }
735*c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
736*c2fc9fa9SBarry Smith #undef __FUNCT__
737*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESReset_VIRS"
738*c2fc9fa9SBarry Smith PetscErrorCode SNESReset_VIRS(SNES snes)
739*c2fc9fa9SBarry Smith {
740*c2fc9fa9SBarry Smith   SNES_VIRS      *vi = (SNES_VIRS*) snes->data;
741*c2fc9fa9SBarry Smith   PetscErrorCode ierr;
742*c2fc9fa9SBarry Smith 
743*c2fc9fa9SBarry Smith   PetscFunctionBegin;
744*c2fc9fa9SBarry Smith   ierr = SNESReset_VI(snes);CHKERRQ(ierr);
745*c2fc9fa9SBarry Smith   ierr = ISDestroy(&vi->IS_inact_prev);CHKERRQ(ierr);
746*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
747*c2fc9fa9SBarry Smith }
748*c2fc9fa9SBarry Smith 
749*c2fc9fa9SBarry Smith 
750*c2fc9fa9SBarry Smith 
751*c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
752*c2fc9fa9SBarry Smith /*MC
753*c2fc9fa9SBarry Smith       SNESVIRS - Reduced space active set solvers for variational inequalities based on Newton's method
754*c2fc9fa9SBarry Smith 
755*c2fc9fa9SBarry Smith 
756*c2fc9fa9SBarry Smith    Level: beginner
757*c2fc9fa9SBarry Smith 
758*c2fc9fa9SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
759*c2fc9fa9SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
760*c2fc9fa9SBarry Smith            SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
761*c2fc9fa9SBarry Smith 
762*c2fc9fa9SBarry Smith M*/
763*c2fc9fa9SBarry Smith EXTERN_C_BEGIN
764*c2fc9fa9SBarry Smith #undef __FUNCT__
765*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESCreate_VIRS"
766*c2fc9fa9SBarry Smith PetscErrorCode  SNESCreate_VIRS(SNES snes)
767*c2fc9fa9SBarry Smith {
768*c2fc9fa9SBarry Smith   PetscErrorCode ierr;
769*c2fc9fa9SBarry Smith   SNES_VIRS      *vi;
770*c2fc9fa9SBarry Smith 
771*c2fc9fa9SBarry Smith   PetscFunctionBegin;
772*c2fc9fa9SBarry Smith   snes->ops->reset           = SNESReset_VIRS;
773*c2fc9fa9SBarry Smith   snes->ops->setup           = SNESSetUp_VIRS;
774*c2fc9fa9SBarry Smith   snes->ops->solve           = SNESSolve_VIRS;
775*c2fc9fa9SBarry Smith   snes->ops->destroy         = SNESDestroy_VI;
776*c2fc9fa9SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_VI;
777*c2fc9fa9SBarry Smith   snes->ops->view            = PETSC_NULL;
778*c2fc9fa9SBarry Smith   snes->ops->converged       = SNESDefaultConverged_VI;
779*c2fc9fa9SBarry Smith 
780*c2fc9fa9SBarry Smith   snes->usesksp             = PETSC_TRUE;
781*c2fc9fa9SBarry Smith   snes->usespc              = PETSC_FALSE;
782*c2fc9fa9SBarry Smith 
783*c2fc9fa9SBarry Smith   ierr                       = PetscNewLog(snes,SNES_VIRS,&vi);CHKERRQ(ierr);
784*c2fc9fa9SBarry Smith   snes->data                 = (void*)vi;
785*c2fc9fa9SBarry Smith   snes->ls_alpha             = 1.e-4;
786*c2fc9fa9SBarry Smith   snes->maxstep              = 1.e8;
787*c2fc9fa9SBarry Smith   snes->steptol              = 1.e-12;
788*c2fc9fa9SBarry Smith   vi->checkredundancy        = PETSC_NULL;
789*c2fc9fa9SBarry Smith 
790*c2fc9fa9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_VI",SNESLineSearchSetType_VI);CHKERRQ(ierr);
791*c2fc9fa9SBarry Smith   ierr = SNESLineSearchSetType(snes, SNES_LS_CUBIC);CHKERRQ(ierr);
792*c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
793*c2fc9fa9SBarry Smith }
794*c2fc9fa9SBarry Smith EXTERN_C_END
795*c2fc9fa9SBarry Smith 
796