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