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