xref: /petsc/src/snes/impls/vi/rs/virs.c (revision e0f5bfbec699682fa3e8b8532b1176849ea4e12a)
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(KSPSetFromOptions(kspnew));*/
291   PetscFunctionReturn(0);
292 }
293 
294 /* Variational Inequality solver using reduce space method. No semismooth algorithm is
295    implemented in this algorithm. It basically identifies the active constraints and does
296    a linear solve on the other variables (those not associated with the active constraints). */
297 PetscErrorCode SNESSolve_VINEWTONRSLS(SNES snes) {
298   SNES_VINEWTONRSLS   *vi = (SNES_VINEWTONRSLS *)snes->data;
299   PetscInt             maxits, i, lits;
300   SNESLineSearchReason lssucceed;
301   PetscReal            fnorm, gnorm, xnorm = 0, ynorm;
302   Vec                  Y, X, F;
303   KSPConvergedReason   kspreason;
304   KSP                  ksp;
305   PC                   pc;
306 
307   PetscFunctionBegin;
308   /* Multigrid must use Galerkin for coarse grids with active set/reduced space methods; cannot rediscretize on coarser grids*/
309   PetscCall(SNESGetKSP(snes, &ksp));
310   PetscCall(KSPGetPC(ksp, &pc));
311   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_BOTH));
312 
313   snes->numFailures            = 0;
314   snes->numLinearSolveFailures = 0;
315   snes->reason                 = SNES_CONVERGED_ITERATING;
316 
317   maxits = snes->max_its;  /* maximum number of iterations */
318   X      = snes->vec_sol;  /* solution vector */
319   F      = snes->vec_func; /* residual vector */
320   Y      = snes->work[0];  /* work vectors */
321 
322   PetscCall(SNESLineSearchSetVIFunctions(snes->linesearch, SNESVIProjectOntoBounds, SNESVIComputeInactiveSetFnorm));
323   PetscCall(SNESLineSearchSetVecs(snes->linesearch, X, NULL, NULL, NULL, NULL));
324   PetscCall(SNESLineSearchSetUp(snes->linesearch));
325 
326   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
327   snes->iter = 0;
328   snes->norm = 0.0;
329   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
330 
331   PetscCall(SNESVIProjectOntoBounds(snes, X));
332   PetscCall(SNESComputeFunction(snes, X, F));
333   PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm));
334   PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- ||x||  */
335   SNESCheckFunctionNorm(snes, fnorm);
336   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
337   snes->norm = fnorm;
338   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
339   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
340   PetscCall(SNESMonitor(snes, 0, fnorm));
341 
342   /* test convergence */
343   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
344   if (snes->reason) PetscFunctionReturn(0);
345 
346   for (i = 0; i < maxits; i++) {
347     IS         IS_act;    /* _act -> active set _inact -> inactive set */
348     IS         IS_redact; /* redundant active set */
349     VecScatter scat_act, scat_inact;
350     PetscInt   nis_act, nis_inact;
351     Vec        Y_act, Y_inact, F_inact;
352     Mat        jac_inact_inact, prejac_inact_inact;
353     PetscBool  isequal;
354 
355     /* Call general purpose update function */
356     PetscTryTypeMethod(snes, update, snes->iter);
357     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
358     SNESCheckJacobianDomainerror(snes);
359 
360     /* Create active and inactive index sets */
361 
362     /*original
363     PetscCall(SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&vi->IS_inact));
364      */
365     PetscCall(SNESVIGetActiveSetIS(snes, X, F, &IS_act));
366 
367     if (vi->checkredundancy) {
368       PetscCall((*vi->checkredundancy)(snes, IS_act, &IS_redact, vi->ctxP));
369       if (IS_redact) {
370         PetscCall(ISSort(IS_redact));
371         PetscCall(ISComplement(IS_redact, X->map->rstart, X->map->rend, &vi->IS_inact));
372         PetscCall(ISDestroy(&IS_redact));
373       } else {
374         PetscCall(ISComplement(IS_act, X->map->rstart, X->map->rend, &vi->IS_inact));
375       }
376     } else {
377       PetscCall(ISComplement(IS_act, X->map->rstart, X->map->rend, &vi->IS_inact));
378     }
379 
380     /* Create inactive set submatrix */
381     PetscCall(MatCreateSubMatrix(snes->jacobian, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &jac_inact_inact));
382 
383     if (0) { /* Dead code (temporary developer hack) */
384       IS keptrows;
385       PetscCall(MatFindNonzeroRows(jac_inact_inact, &keptrows));
386       if (keptrows) {
387         PetscInt        cnt, *nrows, k;
388         const PetscInt *krows, *inact;
389         PetscInt        rstart;
390 
391         PetscCall(MatGetOwnershipRange(jac_inact_inact, &rstart, NULL));
392         PetscCall(MatDestroy(&jac_inact_inact));
393         PetscCall(ISDestroy(&IS_act));
394 
395         PetscCall(ISGetLocalSize(keptrows, &cnt));
396         PetscCall(ISGetIndices(keptrows, &krows));
397         PetscCall(ISGetIndices(vi->IS_inact, &inact));
398         PetscCall(PetscMalloc1(cnt, &nrows));
399         for (k = 0; k < cnt; k++) nrows[k] = inact[krows[k] - rstart];
400         PetscCall(ISRestoreIndices(keptrows, &krows));
401         PetscCall(ISRestoreIndices(vi->IS_inact, &inact));
402         PetscCall(ISDestroy(&keptrows));
403         PetscCall(ISDestroy(&vi->IS_inact));
404 
405         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), cnt, nrows, PETSC_OWN_POINTER, &vi->IS_inact));
406         PetscCall(ISComplement(vi->IS_inact, F->map->rstart, F->map->rend, &IS_act));
407         PetscCall(MatCreateSubMatrix(snes->jacobian, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &jac_inact_inact));
408       }
409     }
410     PetscCall(DMSetVI(snes->dm, vi->IS_inact));
411     /* remove later */
412 
413     /*
414     PetscCall(VecView(vi->xu,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xu))->comm)));
415     PetscCall(VecView(vi->xl,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xl))->comm)));
416     PetscCall(VecView(X,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)X))));
417     PetscCall(VecView(F,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)F))));
418     PetscCall(ISView(vi->IS_inact,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)vi->IS_inact))));
419      */
420 
421     /* Get sizes of active and inactive sets */
422     PetscCall(ISGetLocalSize(IS_act, &nis_act));
423     PetscCall(ISGetLocalSize(vi->IS_inact, &nis_inact));
424 
425     /* Create active and inactive set vectors */
426     PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_inact, &F_inact));
427     PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_act, &Y_act));
428     PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_inact, &Y_inact));
429 
430     /* Create scatter contexts */
431     PetscCall(VecScatterCreate(Y, IS_act, Y_act, NULL, &scat_act));
432     PetscCall(VecScatterCreate(Y, vi->IS_inact, Y_inact, NULL, &scat_inact));
433 
434     /* Do a vec scatter to active and inactive set vectors */
435     PetscCall(VecScatterBegin(scat_inact, F, F_inact, INSERT_VALUES, SCATTER_FORWARD));
436     PetscCall(VecScatterEnd(scat_inact, F, F_inact, INSERT_VALUES, SCATTER_FORWARD));
437 
438     PetscCall(VecScatterBegin(scat_act, Y, Y_act, INSERT_VALUES, SCATTER_FORWARD));
439     PetscCall(VecScatterEnd(scat_act, Y, Y_act, INSERT_VALUES, SCATTER_FORWARD));
440 
441     PetscCall(VecScatterBegin(scat_inact, Y, Y_inact, INSERT_VALUES, SCATTER_FORWARD));
442     PetscCall(VecScatterEnd(scat_inact, Y, Y_inact, INSERT_VALUES, SCATTER_FORWARD));
443 
444     /* Active set direction = 0 */
445     PetscCall(VecSet(Y_act, 0));
446     if (snes->jacobian != snes->jacobian_pre) {
447       PetscCall(MatCreateSubMatrix(snes->jacobian_pre, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &prejac_inact_inact));
448     } else prejac_inact_inact = jac_inact_inact;
449 
450     PetscCall(ISEqual(vi->IS_inact_prev, vi->IS_inact, &isequal));
451     if (!isequal) {
452       PetscCall(SNESVIResetPCandKSP(snes, jac_inact_inact, prejac_inact_inact));
453       PetscCall(PCFieldSplitRestrictIS(pc, vi->IS_inact));
454     }
455 
456     /*      PetscCall(ISView(vi->IS_inact,0)); */
457     /*      PetscCall(ISView(IS_act,0));*/
458     /*      ierr = MatView(snes->jacobian_pre,0); */
459 
460     PetscCall(KSPSetOperators(snes->ksp, jac_inact_inact, prejac_inact_inact));
461     PetscCall(KSPSetUp(snes->ksp));
462     {
463       PC        pc;
464       PetscBool flg;
465       PetscCall(KSPGetPC(snes->ksp, &pc));
466       PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &flg));
467       if (flg) {
468         KSP *subksps;
469         PetscCall(PCFieldSplitGetSubKSP(pc, NULL, &subksps));
470         PetscCall(KSPGetPC(subksps[0], &pc));
471         PetscCall(PetscFree(subksps));
472         PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &flg));
473         if (flg) {
474           PetscInt        n, N = 101 * 101, j, cnts[3] = {0, 0, 0};
475           const PetscInt *ii;
476 
477           PetscCall(ISGetSize(vi->IS_inact, &n));
478           PetscCall(ISGetIndices(vi->IS_inact, &ii));
479           for (j = 0; j < n; j++) {
480             if (ii[j] < N) cnts[0]++;
481             else if (ii[j] < 2 * N) cnts[1]++;
482             else if (ii[j] < 3 * N) cnts[2]++;
483           }
484           PetscCall(ISRestoreIndices(vi->IS_inact, &ii));
485 
486           PetscCall(PCBJacobiSetTotalBlocks(pc, 3, cnts));
487         }
488       }
489     }
490 
491     PetscCall(KSPSolve(snes->ksp, F_inact, Y_inact));
492     PetscCall(VecScatterBegin(scat_act, Y_act, Y, INSERT_VALUES, SCATTER_REVERSE));
493     PetscCall(VecScatterEnd(scat_act, Y_act, Y, INSERT_VALUES, SCATTER_REVERSE));
494     PetscCall(VecScatterBegin(scat_inact, Y_inact, Y, INSERT_VALUES, SCATTER_REVERSE));
495     PetscCall(VecScatterEnd(scat_inact, Y_inact, Y, INSERT_VALUES, SCATTER_REVERSE));
496 
497     PetscCall(VecDestroy(&F_inact));
498     PetscCall(VecDestroy(&Y_act));
499     PetscCall(VecDestroy(&Y_inact));
500     PetscCall(VecScatterDestroy(&scat_act));
501     PetscCall(VecScatterDestroy(&scat_inact));
502     PetscCall(ISDestroy(&IS_act));
503     if (!isequal) {
504       PetscCall(ISDestroy(&vi->IS_inact_prev));
505       PetscCall(ISDuplicate(vi->IS_inact, &vi->IS_inact_prev));
506     }
507     PetscCall(ISDestroy(&vi->IS_inact));
508     PetscCall(MatDestroy(&jac_inact_inact));
509     if (snes->jacobian != snes->jacobian_pre) PetscCall(MatDestroy(&prejac_inact_inact));
510 
511     PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason));
512     if (kspreason < 0) {
513       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
514         PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", number linear solve failures %" PetscInt_FMT " greater than current SNES allowed, stopping solve\n", snes->iter, snes->numLinearSolveFailures));
515         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
516         break;
517       }
518     }
519 
520     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
521     snes->linear_its += lits;
522     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
523     /*
524     if (snes->ops->precheck) {
525       PetscBool changed_y = PETSC_FALSE;
526       PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y);
527     }
528 
529     if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W));
530     */
531     /* Compute a (scaled) negative update in the line search routine:
532          Y <- X - lambda*Y
533        and evaluate G = function(Y) (depends on the line search).
534     */
535     PetscCall(VecCopy(Y, snes->vec_sol_update));
536     ynorm = 1;
537     gnorm = fnorm;
538     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &gnorm, Y));
539     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
540     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm));
541     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lssucceed));
542     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
543     if (snes->domainerror) {
544       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
545       PetscCall(DMDestroyVI(snes->dm));
546       PetscFunctionReturn(0);
547     }
548     if (lssucceed) {
549       if (++snes->numFailures >= snes->maxFailures) {
550         PetscBool ismin;
551         snes->reason = SNES_DIVERGED_LINE_SEARCH;
552         PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, F, X, gnorm, &ismin));
553         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
554         break;
555       }
556     }
557     PetscCall(DMDestroyVI(snes->dm));
558     /* Update function and solution vectors */
559     fnorm = gnorm;
560     /* Monitor convergence */
561     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
562     snes->iter  = i + 1;
563     snes->norm  = fnorm;
564     snes->xnorm = xnorm;
565     snes->ynorm = ynorm;
566     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
567     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
568     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
569     /* Test for convergence, xnorm = || X || */
570     if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
571     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
572     if (snes->reason) break;
573   }
574   /* 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 */
575   PetscCall(DMDestroyVI(snes->dm));
576   if (i == maxits) {
577     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
578     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
579   }
580   PetscFunctionReturn(0);
581 }
582 
583 /*@C
584    SNESVISetRedundancyCheck - Provide a function to check for any redundancy in the VI active set
585 
586    Logically Collective on snes
587 
588    Input Parameters:
589 +  snes - the `SNESVINEWTONRSLS` context
590 .  func - the function to check of redundancies
591 -  ctx - optional context used by the function
592 
593    Level: advanced
594 
595    Note:
596    Sometimes the inactive set will result in a non-singular sub-Jacobian problem that needs to be solved, this allows the user,
597    when they know more about their specific problem to provide a function that removes the redundancy that results in the singular linear system
598 
599 .seealso: `SNESVINEWTONRSLS`, `SNESVIGetInactiveSet()`, `DMSetVI()`
600  @*/
601 PetscErrorCode SNESVISetRedundancyCheck(SNES snes, PetscErrorCode (*func)(SNES, IS, IS *, void *), void *ctx) {
602   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data;
603 
604   PetscFunctionBegin;
605   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
606   vi->checkredundancy = func;
607   vi->ctxP            = ctx;
608   PetscFunctionReturn(0);
609 }
610 
611 #if defined(PETSC_HAVE_MATLAB_ENGINE)
612 #include <engine.h>
613 #include <mex.h>
614 typedef struct {
615   char    *funcname;
616   mxArray *ctx;
617 } SNESMatlabContext;
618 
619 PetscErrorCode SNESVIRedundancyCheck_Matlab(SNES snes, IS is_act, IS *is_redact, void *ctx) {
620   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
621   int                nlhs = 1, nrhs = 5;
622   mxArray           *plhs[1], *prhs[5];
623   long long int      l1 = 0, l2 = 0, ls = 0;
624   PetscInt          *indices = NULL;
625 
626   PetscFunctionBegin;
627   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
628   PetscValidHeaderSpecific(is_act, IS_CLASSID, 2);
629   PetscValidPointer(is_redact, 3);
630   PetscCheckSameComm(snes, 1, is_act, 2);
631 
632   /* Create IS for reduced active set of size 0, its size and indices will
633    bet set by the Matlab function */
634   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), 0, indices, PETSC_OWN_POINTER, is_redact));
635   /* call Matlab function in ctx */
636   PetscCall(PetscArraycpy(&ls, &snes, 1));
637   PetscCall(PetscArraycpy(&l1, &is_act, 1));
638   PetscCall(PetscArraycpy(&l2, is_redact, 1));
639   prhs[0] = mxCreateDoubleScalar((double)ls);
640   prhs[1] = mxCreateDoubleScalar((double)l1);
641   prhs[2] = mxCreateDoubleScalar((double)l2);
642   prhs[3] = mxCreateString(sctx->funcname);
643   prhs[4] = sctx->ctx;
644   PetscCall(mexCallMATLAB(nlhs, plhs, nrhs, prhs, "PetscSNESVIRedundancyCheckInternal"));
645   PetscCall(mxGetScalar(plhs[0]));
646   mxDestroyArray(prhs[0]);
647   mxDestroyArray(prhs[1]);
648   mxDestroyArray(prhs[2]);
649   mxDestroyArray(prhs[3]);
650   mxDestroyArray(plhs[0]);
651   PetscFunctionReturn(0);
652 }
653 
654 PetscErrorCode SNESVISetRedundancyCheckMatlab(SNES snes, const char *func, mxArray *ctx) {
655   SNESMatlabContext *sctx;
656 
657   PetscFunctionBegin;
658   /* currently sctx is memory bleed */
659   PetscCall(PetscNew(&sctx));
660   PetscCall(PetscStrallocpy(func, &sctx->funcname));
661   sctx->ctx = mxDuplicateArray(ctx);
662   PetscCall(SNESVISetRedundancyCheck(snes, SNESVIRedundancyCheck_Matlab, sctx));
663   PetscFunctionReturn(0);
664 }
665 
666 #endif
667 
668 /*
669    SNESSetUp_VINEWTONRSLS - Sets up the internal data structures for the later use
670    of the SNESVI nonlinear solver.
671 
672    Input Parameter:
673 .  snes - the SNES context
674 
675    Application Interface Routine: SNESSetUp()
676 
677    Note:
678    For basic use of the SNES solvers, the user need not explicitly call
679    SNESSetUp(), since these actions will automatically occur during
680    the call to SNESSolve().
681  */
682 PetscErrorCode SNESSetUp_VINEWTONRSLS(SNES snes) {
683   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data;
684   PetscInt          *indices;
685   PetscInt           i, n, rstart, rend;
686   SNESLineSearch     linesearch;
687 
688   PetscFunctionBegin;
689   PetscCall(SNESSetUp_VI(snes));
690 
691   /* Set up previous active index set for the first snes solve
692    vi->IS_inact_prev = 0,1,2,....N */
693 
694   PetscCall(VecGetOwnershipRange(snes->vec_sol, &rstart, &rend));
695   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
696   PetscCall(PetscMalloc1(n, &indices));
697   for (i = 0; i < n; i++) indices[i] = rstart + i;
698   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), n, indices, PETSC_OWN_POINTER, &vi->IS_inact_prev));
699 
700   /* set the line search functions */
701   if (!snes->linesearch) {
702     PetscCall(SNESGetLineSearch(snes, &linesearch));
703     PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
704   }
705   PetscFunctionReturn(0);
706 }
707 
708 PetscErrorCode SNESReset_VINEWTONRSLS(SNES snes) {
709   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data;
710 
711   PetscFunctionBegin;
712   PetscCall(SNESReset_VI(snes));
713   PetscCall(ISDestroy(&vi->IS_inact_prev));
714   PetscFunctionReturn(0);
715 }
716 
717 /*MC
718       SNESVINEWTONRSLS - Reduced space active set solvers for variational inequalities based on Newton's method
719 
720    Options Database Keys:
721 +   -snes_type <vinewtonssls,vinewtonrsls> - a semi-smooth solver, a reduced space active set method
722 -   -snes_vi_monitor - prints the number of active constraints at each iteration.
723 
724    Level: beginner
725 
726    References:
727 .  * - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale
728      Applications, Optimization Methods and Software, 21 (2006).
729 
730    Note:
731    At each set of this methods the algorithm produces an inactive set of variables that are constrained to their current values
732    (because changing these values would result in those variables no longer satisfying the inequality constraints)
733    and produces a step direction by solving the linear system arising from the Jacobian with the inactive variables removed. In other
734    words on a reduced space of the solution space. Based on the Newton update it then adjusts the inactive sep for the next iteration.
735 
736 .seealso: `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONSSLS`, `SNESNEWTONTRDC`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESVIGetInactiveSet()`, `DMSetVI()`, `SNESVISetRedundancyCheck()`
737 M*/
738 PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONRSLS(SNES snes) {
739   SNES_VINEWTONRSLS *vi;
740   SNESLineSearch     linesearch;
741 
742   PetscFunctionBegin;
743   snes->ops->reset          = SNESReset_VINEWTONRSLS;
744   snes->ops->setup          = SNESSetUp_VINEWTONRSLS;
745   snes->ops->solve          = SNESSolve_VINEWTONRSLS;
746   snes->ops->destroy        = SNESDestroy_VI;
747   snes->ops->setfromoptions = SNESSetFromOptions_VI;
748   snes->ops->view           = NULL;
749   snes->ops->converged      = SNESConvergedDefault_VI;
750 
751   snes->usesksp = PETSC_TRUE;
752   snes->usesnpc = PETSC_FALSE;
753 
754   PetscCall(SNESGetLineSearch(snes, &linesearch));
755   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
756   PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0));
757 
758   snes->alwayscomputesfinalresidual = PETSC_TRUE;
759 
760   PetscCall(PetscNew(&vi));
761   snes->data          = (void *)vi;
762   vi->checkredundancy = NULL;
763 
764   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI));
765   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI));
766   PetscFunctionReturn(0);
767 }
768