Actual source code: virs.c

  1: #include <../src/snes/impls/vi/rs/virsimpl.h>
  2: #include <petsc/private/dmimpl.h>
  3: #include <petsc/private/vecimpl.h>

  5: /*@
  6:   SNESVIGetInactiveSet - Gets the global indices for the inactive set variables (these correspond to the degrees of freedom the linear
  7:   system is solved on)

  9:   Input Parameter:
 10: . snes - the `SNES` context

 12:   Output Parameter:
 13: . inact - inactive set index set

 15:   Level: advanced

 17:   Note:
 18:   See `SNESVINEWTONRSLS` for a concise description of the active and inactive sets

 20: .seealso: [](ch_snes), `SNES`, `SNESVINEWTONRSLS`
 21: @*/
 22: PetscErrorCode SNESVIGetInactiveSet(SNES snes, IS *inact)
 23: {
 24:   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data;

 26:   PetscFunctionBegin;
 27:   *inact = vi->IS_inact;
 28:   PetscFunctionReturn(PETSC_SUCCESS);
 29: }

 31: /*
 32:     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
 33:   defined by the reduced space method.

 35:     Simple calls the regular DM interpolation and restricts it to operation on the variables not associated with active constraints.
 36: */
 37: typedef struct {
 38:   PetscInt n; /* size of vectors in the reduced DM space */
 39:   IS       inactive;

 41:   PetscErrorCode (*createinterpolation)(DM, DM, Mat *, Vec *); /* DM's original routines */
 42:   PetscErrorCode (*coarsen)(DM, MPI_Comm, DM *);
 43:   PetscErrorCode (*createglobalvector)(DM, Vec *);
 44:   PetscErrorCode (*createinjection)(DM, DM, Mat *);
 45:   PetscErrorCode (*hascreateinjection)(DM, PetscBool *);

 47:   DM dm; /* when destroying this object we need to reset the above function into the base DM */
 48: } DM_SNESVI;

 50: /*
 51:      DMCreateGlobalVector_SNESVI - Creates global vector of the size of the reduced space
 52: */
 53: static PetscErrorCode DMCreateGlobalVector_SNESVI(DM dm, Vec *vec)
 54: {
 55:   PetscContainer isnes;
 56:   DM_SNESVI     *dmsnesvi;

 58:   PetscFunctionBegin;
 59:   PetscCall(PetscObjectQuery((PetscObject)dm, "VI", (PetscObject *)&isnes));
 60:   PetscCheck(isnes, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Composed SNES is missing");
 61:   PetscCall(PetscContainerGetPointer(isnes, &dmsnesvi));
 62:   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)dm), dmsnesvi->n, PETSC_DETERMINE, vec));
 63:   PetscCall(VecSetDM(*vec, dm));
 64:   PetscFunctionReturn(PETSC_SUCCESS);
 65: }

 67: /*
 68:      DMCreateInterpolation_SNESVI - Modifieds the interpolation obtained from the DM by removing all rows and columns associated with active constraints.
 69: */
 70: static PetscErrorCode DMCreateInterpolation_SNESVI(DM dm1, DM dm2, Mat *mat, Vec *vec)
 71: {
 72:   PetscContainer isnes;
 73:   DM_SNESVI     *dmsnesvi1, *dmsnesvi2;
 74:   Mat            interp;

 76:   PetscFunctionBegin;
 77:   PetscCall(PetscObjectQuery((PetscObject)dm1, "VI", (PetscObject *)&isnes));
 78:   PetscCheck(isnes, PetscObjectComm((PetscObject)dm1), PETSC_ERR_PLIB, "Composed VI data structure is missing");
 79:   PetscCall(PetscContainerGetPointer(isnes, &dmsnesvi1));
 80:   PetscCall(PetscObjectQuery((PetscObject)dm2, "VI", (PetscObject *)&isnes));
 81:   PetscCheck(isnes, PetscObjectComm((PetscObject)dm2), PETSC_ERR_PLIB, "Composed VI data structure is missing");
 82:   PetscCall(PetscContainerGetPointer(isnes, &dmsnesvi2));

 84:   PetscCall((*dmsnesvi1->createinterpolation)(dm1, dm2, &interp, NULL));
 85:   PetscCall(MatCreateSubMatrix(interp, dmsnesvi2->inactive, dmsnesvi1->inactive, MAT_INITIAL_MATRIX, mat));
 86:   PetscCall(MatDestroy(&interp));
 87:   *vec = NULL;
 88:   PetscFunctionReturn(PETSC_SUCCESS);
 89: }

 91: /*
 92:      DMCoarsen_SNESVI - Computes the regular coarsened DM then computes additional information about its inactive set
 93: */
 94: static PetscErrorCode DMCoarsen_SNESVI(DM dm1, MPI_Comm comm, DM *dm2)
 95: {
 96:   PetscContainer  isnes;
 97:   DM_SNESVI      *dmsnesvi1;
 98:   Vec             finemarked, coarsemarked;
 99:   IS              inactive;
100:   Mat             inject;
101:   const PetscInt *index;
102:   PetscInt        n, k, cnt = 0, rstart, *coarseindex;
103:   PetscScalar    *marked;

105:   PetscFunctionBegin;
106:   PetscCall(PetscObjectQuery((PetscObject)dm1, "VI", (PetscObject *)&isnes));
107:   PetscCheck(isnes, PetscObjectComm((PetscObject)dm1), PETSC_ERR_PLIB, "Composed VI data structure is missing");
108:   PetscCall(PetscContainerGetPointer(isnes, &dmsnesvi1));

110:   /* get the original coarsen */
111:   PetscCall((*dmsnesvi1->coarsen)(dm1, comm, dm2));

113:   /* not sure why this extra reference is needed, but without the dm2 disappears too early */
114:   /* Updating the KSPCreateVecs() to avoid using DMGetGlobalVector() when matrix is available removes the need for this reference? */
115:   /*  PetscCall(PetscObjectReference((PetscObject)*dm2));*/

117:   /* need to set back global vectors in order to use the original injection */
118:   PetscCall(DMClearGlobalVectors(dm1));

120:   dm1->ops->createglobalvector = dmsnesvi1->createglobalvector;

122:   PetscCall(DMCreateGlobalVector(dm1, &finemarked));
123:   PetscCall(DMCreateGlobalVector(*dm2, &coarsemarked));

125:   /*
126:      fill finemarked with locations of inactive points
127:   */
128:   PetscCall(ISGetIndices(dmsnesvi1->inactive, &index));
129:   PetscCall(ISGetLocalSize(dmsnesvi1->inactive, &n));
130:   for (k = 0; k < n; k++) PetscCall(VecSetValue(finemarked, index[k], 1.0, INSERT_VALUES));
131:   PetscCall(VecAssemblyBegin(finemarked));
132:   PetscCall(VecAssemblyEnd(finemarked));

134:   PetscCall(DMCreateInjection(*dm2, dm1, &inject));
135:   PetscCall(MatRestrict(inject, finemarked, coarsemarked));
136:   PetscCall(MatDestroy(&inject));

138:   /*
139:      create index set list of coarse inactive points from coarsemarked
140:   */
141:   PetscCall(VecGetLocalSize(coarsemarked, &n));
142:   PetscCall(VecGetOwnershipRange(coarsemarked, &rstart, NULL));
143:   PetscCall(VecGetArray(coarsemarked, &marked));
144:   for (k = 0; k < n; k++) {
145:     if (marked[k] != 0.0) cnt++;
146:   }
147:   PetscCall(PetscMalloc1(cnt, &coarseindex));
148:   cnt = 0;
149:   for (k = 0; k < n; k++) {
150:     if (marked[k] != 0.0) coarseindex[cnt++] = k + rstart;
151:   }
152:   PetscCall(VecRestoreArray(coarsemarked, &marked));
153:   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)coarsemarked), cnt, coarseindex, PETSC_OWN_POINTER, &inactive));

155:   PetscCall(DMClearGlobalVectors(dm1));

157:   dm1->ops->createglobalvector = DMCreateGlobalVector_SNESVI;

159:   PetscCall(DMSetVI(*dm2, inactive));

161:   PetscCall(VecDestroy(&finemarked));
162:   PetscCall(VecDestroy(&coarsemarked));
163:   PetscCall(ISDestroy(&inactive));
164:   PetscFunctionReturn(PETSC_SUCCESS);
165: }

167: static PetscErrorCode DMDestroy_SNESVI(PetscCtxRt ctx)
168: {
169:   DM_SNESVI *dmsnesvi = *(DM_SNESVI **)ctx;

171:   PetscFunctionBegin;
172:   /* reset the base methods in the DM object that were changed when the DM_SNESVI was reset */
173:   dmsnesvi->dm->ops->createinterpolation = dmsnesvi->createinterpolation;
174:   dmsnesvi->dm->ops->coarsen             = dmsnesvi->coarsen;
175:   dmsnesvi->dm->ops->createglobalvector  = dmsnesvi->createglobalvector;
176:   dmsnesvi->dm->ops->createinjection     = dmsnesvi->createinjection;
177:   dmsnesvi->dm->ops->hascreateinjection  = dmsnesvi->hascreateinjection;
178:   /* need to clear out this vectors because some of them may not have a reference to the DM
179:     but they are counted as having references to the DM in DMDestroy() */
180:   PetscCall(DMClearGlobalVectors(dmsnesvi->dm));

182:   PetscCall(ISDestroy(&dmsnesvi->inactive));
183:   PetscCall(PetscFree(dmsnesvi));
184:   PetscFunctionReturn(PETSC_SUCCESS);
185: }

187: /*@
188:   DMSetVI - Marks a `DM` as associated with a VI problem. This causes the interpolation/restriction operators to
189:   be restricted to only those variables NOT associated with active constraints.

191:   Logically Collective

193:   Input Parameters:
194: + dm       - the `DM` object
195: - inactive - an `IS` indicating which points are currently not active

197:   Level: intermediate

199: .seealso: [](ch_snes), `SNES`, `SNESVINEWTONRSLS`, `SNESVIGetInactiveSet()`
200: @*/
201: PetscErrorCode DMSetVI(DM dm, IS inactive)
202: {
203:   PetscContainer isnes;
204:   DM_SNESVI     *dmsnesvi;

206:   PetscFunctionBegin;
207:   if (!dm) PetscFunctionReturn(PETSC_SUCCESS);

209:   PetscCall(PetscObjectReference((PetscObject)inactive));

211:   PetscCall(PetscObjectQuery((PetscObject)dm, "VI", (PetscObject *)&isnes));
212:   if (!isnes) {
213:     PetscCall(PetscNew(&dmsnesvi));
214:     PetscCall(PetscObjectContainerCompose((PetscObject)dm, "VI", dmsnesvi, DMDestroy_SNESVI));

216:     dmsnesvi->createinterpolation = dm->ops->createinterpolation;
217:     dm->ops->createinterpolation  = DMCreateInterpolation_SNESVI;
218:     dmsnesvi->coarsen             = dm->ops->coarsen;
219:     dm->ops->coarsen              = DMCoarsen_SNESVI;
220:     dmsnesvi->createglobalvector  = dm->ops->createglobalvector;
221:     dm->ops->createglobalvector   = DMCreateGlobalVector_SNESVI;
222:     dmsnesvi->createinjection     = dm->ops->createinjection;
223:     dm->ops->createinjection      = NULL;
224:     dmsnesvi->hascreateinjection  = dm->ops->hascreateinjection;
225:     dm->ops->hascreateinjection   = NULL;
226:   } else {
227:     PetscCall(PetscContainerGetPointer(isnes, &dmsnesvi));
228:     PetscCall(ISDestroy(&dmsnesvi->inactive));
229:   }
230:   PetscCall(DMClearGlobalVectors(dm));
231:   PetscCall(ISGetLocalSize(inactive, &dmsnesvi->n));

233:   dmsnesvi->inactive = inactive;
234:   dmsnesvi->dm       = dm;
235:   PetscFunctionReturn(PETSC_SUCCESS);
236: }

238: /*
239:      DMDestroyVI - Frees the DM_SNESVI object contained in the DM
240:          - also resets the function pointers in the DM for createinterpolation() etc to use the original DM
241: */
242: PetscErrorCode DMDestroyVI(DM dm)
243: {
244:   PetscFunctionBegin;
245:   if (!dm) PetscFunctionReturn(PETSC_SUCCESS);
246:   PetscCall(PetscObjectCompose((PetscObject)dm, "VI", (PetscObject)NULL));
247:   PetscFunctionReturn(PETSC_SUCCESS);
248: }

250: /* Create active and inactive set vectors. The local size of this vector is set and PETSc computes the global size */
251: static PetscErrorCode SNESCreateSubVectors_VINEWTONRSLS(SNES snes, PetscInt n, Vec *newv)
252: {
253:   Vec v;

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(PETSC_SUCCESS);
261: }

263: /* Resets the snes PC and KSP when the active set sizes change */
264: static PetscErrorCode SNESVIResetPCandKSP(SNES snes, Mat Amat, Mat Pmat)
265: {
266:   KSP snesksp;

268:   PetscFunctionBegin;
269:   PetscCall(SNESGetKSP(snes, &snesksp));
270:   PetscCall(KSPReset(snesksp));
271:   PetscCall(KSPResetFromOptions(snesksp));

273:   /*
274:   KSP                    kspnew;
275:   PC                     pcnew;
276:   MatSolverType          stype;

278:   PetscCall(KSPCreate(PetscObjectComm((PetscObject)snes),&kspnew));
279:   kspnew->pc_side = snesksp->pc_side;
280:   kspnew->rtol    = snesksp->rtol;
281:   kspnew->abstol    = snesksp->abstol;
282:   kspnew->max_it  = snesksp->max_it;
283:   PetscCall(KSPSetType(kspnew,((PetscObject)snesksp)->type_name));
284:   PetscCall(KSPGetPC(kspnew,&pcnew));
285:   PetscCall(PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name));
286:   PetscCall(PCSetOperators(kspnew->pc,Amat,Pmat));
287:   PetscCall(PCFactorGetMatSolverType(snesksp->pc,&stype));
288:   PetscCall(PCFactorSetMatSolverType(kspnew->pc,stype));
289:   PetscCall(KSPDestroy(&snesksp));
290:   snes->ksp = kspnew;
291:    PetscCall(KSPSetFromOptions(kspnew));*/
292:   PetscFunctionReturn(PETSC_SUCCESS);
293: }

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: static PetscErrorCode SNESSolve_VINEWTONRSLS(SNES snes)
299: {
300:   SNES_VINEWTONRSLS   *vi = (SNES_VINEWTONRSLS *)snes->data;
301:   PetscInt             maxits, i, lits;
302:   SNESLineSearchReason lsreason;
303:   PetscReal            fnorm, gnorm, xnorm = 0, ynorm;
304:   Vec                  Y, X, F;
305:   KSPConvergedReason   kspreason;
306:   KSP                  ksp;
307:   PC                   pc;

309:   PetscFunctionBegin;
310:   /* Multigrid must use Galerkin for coarse grids with active set/reduced space methods; cannot rediscretize on coarser grids*/
311:   PetscCall(SNESGetKSP(snes, &ksp));
312:   PetscCall(KSPGetPC(ksp, &pc));
313:   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_BOTH));

315:   snes->numFailures            = 0;
316:   snes->numLinearSolveFailures = 0;
317:   snes->reason                 = SNES_CONVERGED_ITERATING;

319:   maxits = snes->max_its;  /* maximum number of iterations */
320:   X      = snes->vec_sol;  /* solution vector */
321:   F      = snes->vec_func; /* residual vector */
322:   Y      = snes->work[0];  /* work vectors */

324:   PetscCall(SNESLineSearchSetVIFunctions(snes->linesearch, SNESVIProjectOntoBounds, SNESVIComputeInactiveSetFnorm, SNESVIComputeInactiveSetFtY));
325:   PetscCall(SNESLineSearchSetVecs(snes->linesearch, X, NULL, NULL, NULL, NULL));
326:   PetscCall(SNESLineSearchSetUp(snes->linesearch));

328:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
329:   snes->iter = 0;
330:   snes->norm = 0.0;
331:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));

333:   PetscCall(SNESVIProjectOntoBounds(snes, X));
334:   PetscCall(SNESComputeFunction(snes, X, F));
335:   PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm));
336:   PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- ||x||  */
337:   SNESCheckFunctionDomainError(snes, fnorm);
338:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
339:   snes->norm = fnorm;
340:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
341:   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));

343:   /* test convergence */
344:   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
345:   PetscCall(SNESMonitor(snes, 0, fnorm));
346:   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);

348:   for (i = 0; i < maxits; i++) {
349:     IS         IS_act;    /* _act -> active set _inact -> inactive set */
350:     IS         IS_redact; /* redundant active set */
351:     VecScatter scat_act, scat_inact;
352:     PetscInt   nis_act, nis_inact;
353:     Vec        Y_act, Y_inact, F_inact;
354:     Mat        jac_inact_inact, prejac_inact_inact;
355:     PetscBool  isequal;

357:     /* Call general purpose update function */
358:     PetscTryTypeMethod(snes, update, snes->iter);
359:     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
360:     SNESCheckJacobianDomainError(snes);

362:     /* Create active and inactive index sets */

364:     /*original
365:     PetscCall(SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&vi->IS_inact));
366:      */
367:     PetscCall(SNESVIGetActiveSetIS(snes, X, F, &IS_act));

369:     if (vi->checkredundancy) {
370:       PetscCall((*vi->checkredundancy)(snes, IS_act, &IS_redact, vi->ctxP));
371:       if (IS_redact) {
372:         PetscCall(ISSort(IS_redact));
373:         PetscCall(ISComplement(IS_redact, X->map->rstart, X->map->rend, &vi->IS_inact));
374:         PetscCall(ISDestroy(&IS_redact));
375:       } else {
376:         PetscCall(ISComplement(IS_act, X->map->rstart, X->map->rend, &vi->IS_inact));
377:       }
378:     } else {
379:       PetscCall(ISComplement(IS_act, X->map->rstart, X->map->rend, &vi->IS_inact));
380:     }

382:     /* Create inactive set submatrix */
383:     PetscCall(MatCreateSubMatrix(snes->jacobian, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &jac_inact_inact));

385:     if (0) { /* Dead code (temporary developer hack) */
386:       IS keptrows;
387:       PetscCall(MatFindNonzeroRows(jac_inact_inact, &keptrows));
388:       if (keptrows) {
389:         PetscInt        cnt, *nrows, k;
390:         const PetscInt *krows, *inact;
391:         PetscInt        rstart;

393:         PetscCall(MatGetOwnershipRange(jac_inact_inact, &rstart, NULL));
394:         PetscCall(MatDestroy(&jac_inact_inact));
395:         PetscCall(ISDestroy(&IS_act));

397:         PetscCall(ISGetLocalSize(keptrows, &cnt));
398:         PetscCall(ISGetIndices(keptrows, &krows));
399:         PetscCall(ISGetIndices(vi->IS_inact, &inact));
400:         PetscCall(PetscMalloc1(cnt, &nrows));
401:         for (k = 0; k < cnt; k++) nrows[k] = inact[krows[k] - rstart];
402:         PetscCall(ISRestoreIndices(keptrows, &krows));
403:         PetscCall(ISRestoreIndices(vi->IS_inact, &inact));
404:         PetscCall(ISDestroy(&keptrows));
405:         PetscCall(ISDestroy(&vi->IS_inact));

407:         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), cnt, nrows, PETSC_OWN_POINTER, &vi->IS_inact));
408:         PetscCall(ISComplement(vi->IS_inact, F->map->rstart, F->map->rend, &IS_act));
409:         PetscCall(MatCreateSubMatrix(snes->jacobian, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &jac_inact_inact));
410:       }
411:     }
412:     PetscCall(DMSetVI(snes->dm, vi->IS_inact));
413:     /* remove later */

415:     /*
416:     PetscCall(VecView(vi->xu,PETSC_VIEWER_BINARY_(((PetscObject)vi->xu)->comm)));
417:     PetscCall(VecView(vi->xl,PETSC_VIEWER_BINARY_(((PetscObject)vi->xl)->comm)));
418:     PetscCall(VecView(X,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)X))));
419:     PetscCall(VecView(F,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)F))));
420:     PetscCall(ISView(vi->IS_inact,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)vi->IS_inact))));
421:      */

423:     /* Get sizes of active and inactive sets */
424:     PetscCall(ISGetLocalSize(IS_act, &nis_act));
425:     PetscCall(ISGetLocalSize(vi->IS_inact, &nis_inact));

427:     /* Create active and inactive set vectors */
428:     PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_inact, &F_inact));
429:     PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_act, &Y_act));
430:     PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_inact, &Y_inact));

432:     /* Create scatter contexts */
433:     PetscCall(VecScatterCreate(Y, IS_act, Y_act, NULL, &scat_act));
434:     PetscCall(VecScatterCreate(Y, vi->IS_inact, Y_inact, NULL, &scat_inact));

436:     /* Do a vec scatter to active and inactive set vectors */
437:     PetscCall(VecScatterBegin(scat_inact, F, F_inact, INSERT_VALUES, SCATTER_FORWARD));
438:     PetscCall(VecScatterEnd(scat_inact, F, F_inact, INSERT_VALUES, SCATTER_FORWARD));

440:     PetscCall(VecScatterBegin(scat_act, Y, Y_act, INSERT_VALUES, SCATTER_FORWARD));
441:     PetscCall(VecScatterEnd(scat_act, Y, Y_act, INSERT_VALUES, SCATTER_FORWARD));

443:     PetscCall(VecScatterBegin(scat_inact, Y, Y_inact, INSERT_VALUES, SCATTER_FORWARD));
444:     PetscCall(VecScatterEnd(scat_inact, Y, Y_inact, INSERT_VALUES, SCATTER_FORWARD));

446:     /* Active set direction = 0 */
447:     PetscCall(VecSet(Y_act, 0));
448:     if (snes->jacobian != snes->jacobian_pre) 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;

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:     }

457:     /*      PetscCall(ISView(vi->IS_inact,0)); */
458:     /*      PetscCall(ISView(IS_act,0));*/
459:     /*      ierr = MatView(snes->jacobian_pre,0); */

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;

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));

487:           PetscCall(PCBJacobiSetTotalBlocks(pc, 3, cnts));
488:         }
489:       }
490:     }

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));

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));

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:     }

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:     }

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(DMDestroyVI(snes->dm));
541:     if (snes->reason) break;
542:     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsreason));
543:     if (lsreason) {
544:       if (snes->stol * xnorm > ynorm) {
545:         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
546:         break;
547:       } else if (lsreason == SNES_LINESEARCH_FAILED_FUNCTION_DOMAIN) {
548:         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
549:         break;
550:       } else if (lsreason == SNES_LINESEARCH_FAILED_NANORINF) {
551:         snes->reason = SNES_DIVERGED_FUNCTION_NANORINF;
552:         break;
553:       } else if (lsreason == SNES_LINESEARCH_FAILED_OBJECTIVE_DOMAIN) {
554:         snes->reason = SNES_DIVERGED_OBJECTIVE_DOMAIN;
555:         break;
556:       } else if (lsreason == SNES_LINESEARCH_FAILED_JACOBIAN_DOMAIN) {
557:         snes->reason = SNES_DIVERGED_JACOBIAN_DOMAIN;
558:         break;
559:       } else if (++snes->numFailures >= snes->maxFailures) {
560:         PetscBool ismin;

562:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
563:         PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, F, X, gnorm, &ismin));
564:         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
565:         break;
566:       }
567:     }
568:     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm));
569:     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lsreason));
570:     /* Update function and solution vectors */
571:     fnorm = gnorm;
572:     /* Monitor convergence */
573:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
574:     snes->iter  = i + 1;
575:     snes->norm  = fnorm;
576:     snes->xnorm = xnorm;
577:     snes->ynorm = ynorm;
578:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
579:     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
580:     /* Test for convergence, xnorm = || X || */
581:     if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
582:     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
583:     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
584:     if (snes->reason) break;
585:   }
586:   /* 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 */
587:   PetscCall(DMDestroyVI(snes->dm));
588:   PetscFunctionReturn(PETSC_SUCCESS);
589: }

591: /*@C
592:   SNESVISetRedundancyCheck - Provide a function to check for any redundancy in the VI active set

594:   Logically Collective

596:   Input Parameters:
597: + snes - the `SNESVINEWTONRSLS` context
598: . func - the function to check of redundancies
599: - ctx  - optional context used by the function

601:   Calling sequence of func:
602: + snes      - the `SNES` context
603: . is_act    - the set of points in the active sets
604: . is_redact - output, the set of points in the non-redundant active set
605: - ctx       - optional context

607:   Level: advanced

609:   Note:
610:   Sometimes the inactive set will result in a singular sub-Jacobian problem that needs to be solved, this allows the user,
611:   when they know more about their specific problem to provide a function that removes the redundancy that results in the singular linear system

613:   See `SNESVINEWTONRSLS` for a concise description of the active and inactive sets

615: .seealso: [](ch_snes), `SNES`, `SNESVINEWTONRSLS`, `SNESVIGetInactiveSet()`, `DMSetVI()`
616:  @*/
617: PetscErrorCode SNESVISetRedundancyCheck(SNES snes, PetscErrorCode (*func)(SNES snes, IS is_act, IS *is_redact, PetscCtx ctx), PetscCtx ctx)
618: {
619:   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data;

621:   PetscFunctionBegin;
623:   vi->checkredundancy = func;
624:   vi->ctxP            = ctx;
625:   PetscFunctionReturn(PETSC_SUCCESS);
626: }

628: #if defined(PETSC_HAVE_MATLAB)
629:   #include <engine.h>
630:   #include <mex.h>
631: typedef struct {
632:   char    *funcname;
633:   mxArray *ctx;
634: } SNESMatlabContext;

636: PetscErrorCode SNESVIRedundancyCheck_Matlab(SNES snes, IS is_act, IS *is_redact, PetscCtx ctx)
637: {
638:   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
639:   int                nlhs = 1, nrhs = 5;
640:   mxArray           *plhs[1], *prhs[5];
641:   long long int      l1 = 0, l2 = 0, ls = 0;
642:   PetscInt          *indices = NULL;

644:   PetscFunctionBegin;
647:   PetscAssertPointer(is_redact, 3);
648:   PetscCheckSameComm(snes, 1, is_act, 2);

650:   /* Create IS for reduced active set of size 0, its size and indices will
651:    bet set by the MATLAB function */
652:   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), 0, indices, PETSC_OWN_POINTER, is_redact));
653:   /* call MATLAB function in ctx */
654:   PetscCall(PetscArraycpy(&ls, &snes, 1));
655:   PetscCall(PetscArraycpy(&l1, &is_act, 1));
656:   PetscCall(PetscArraycpy(&l2, is_redact, 1));
657:   prhs[0] = mxCreateDoubleScalar((double)ls);
658:   prhs[1] = mxCreateDoubleScalar((double)l1);
659:   prhs[2] = mxCreateDoubleScalar((double)l2);
660:   prhs[3] = mxCreateString(sctx->funcname);
661:   prhs[4] = sctx->ctx;
662:   PetscCall(mexCallMATLAB(nlhs, plhs, nrhs, prhs, "PetscSNESVIRedundancyCheckInternal"));
663:   PetscCall(mxGetScalar(plhs[0]));
664:   mxDestroyArray(prhs[0]);
665:   mxDestroyArray(prhs[1]);
666:   mxDestroyArray(prhs[2]);
667:   mxDestroyArray(prhs[3]);
668:   mxDestroyArray(plhs[0]);
669:   PetscFunctionReturn(PETSC_SUCCESS);
670: }

672: PetscErrorCode SNESVISetRedundancyCheckMatlab(SNES snes, const char *func, mxArray *ctx)
673: {
674:   SNESMatlabContext *sctx;

676:   PetscFunctionBegin;
677:   /* currently sctx is memory bleed */
678:   PetscCall(PetscNew(&sctx));
679:   PetscCall(PetscStrallocpy(func, &sctx->funcname));
680:   sctx->ctx = mxDuplicateArray(ctx);
681:   PetscCall(SNESVISetRedundancyCheck(snes, SNESVIRedundancyCheck_Matlab, sctx));
682:   PetscFunctionReturn(PETSC_SUCCESS);
683: }

685: #endif

687: static PetscErrorCode SNESSetUp_VINEWTONRSLS(SNES snes)
688: {
689:   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data;
690:   PetscInt          *indices;
691:   PetscInt           i, n, rstart, rend;
692:   SNESLineSearch     linesearch;

694:   PetscFunctionBegin;
695:   PetscCall(SNESSetUp_VI(snes));

697:   /* Set up previous active index set for the first snes solve
698:    vi->IS_inact_prev = 0,1,2,....N */

700:   PetscCall(VecGetOwnershipRange(snes->work[0], &rstart, &rend));
701:   PetscCall(VecGetLocalSize(snes->work[0], &n));
702:   PetscCall(PetscMalloc1(n, &indices));
703:   for (i = 0; i < n; i++) indices[i] = rstart + i;
704:   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), n, indices, PETSC_OWN_POINTER, &vi->IS_inact_prev));

706:   /* set the line search functions */
707:   if (!snes->linesearch) {
708:     PetscCall(SNESGetLineSearch(snes, &linesearch));
709:     PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
710:   }
711:   PetscFunctionReturn(PETSC_SUCCESS);
712: }

714: static PetscErrorCode SNESReset_VINEWTONRSLS(SNES snes)
715: {
716:   SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data;

718:   PetscFunctionBegin;
719:   PetscCall(SNESReset_VI(snes));
720:   PetscCall(ISDestroy(&vi->IS_inact_prev));
721:   PetscFunctionReturn(PETSC_SUCCESS);
722: }

724: /*MC
725:    SNESVINEWTONRSLS - Reduced space active set solvers for variational inequalities based on Newton's method

727:    Options Database Keys:
728: +  -snes_type (vinewtonssls|vinewtonrsls) - A semi-smooth solver or a reduced space active set method
729: .  -snes_vi_zero_tolerance                - Tolerance for considering $u_i$ value to be on a bound.
730: .  -snes_vi_monitor                       - Prints the number of active constraints (inactive set points) at each iteration.
731: .  -snes_vi_monitor_residual              - View the residual vector at each iteration, using zero for active constraints (i.e. the inactive variables).
732: -  -snes_vi_monitor_active                - View the active set by outputting a one for vector components in the active set and zero for the inactive.

734:    Level: beginner

736:    Note:
737:    Reduced-space (active set methods) work as follows at each iteration\:
738:    - The algorithm produces an inactive set of variables, that is a list of variables whose values will not be changed in the current iteration, i.e. they
739:      are to be constrained to their current values. These are all the variables that are on the lower bound, that is $u_i = L_i$ with also $[F(u)]_i \ge 0$ or
740:      the upper bound $u_i = U_i$ with also $[F(u)]_i \le 0.$
741:    - A step direction is obtained by solving the linear system arising from the Jacobian used in Newton's method but with the inactive variables removed
742:      from both the rows and columns.
743:    - A line search is then used to update the active variables (the inactive set of variables are not changed).

745:    The inactive set is chosen based on the sign of $[F(u)]_i$ because this gives exactly the set of points that would be moved outside of the domain given
746:    an infinitesimal Newton (or even Richardson) step and our goal is to remain within the bounds, that is, to continue to satisfy the inequality constraints.

748:    See {cite}`benson2006flexible`

750: .seealso: [](ch_snes), `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONSSLS`, `SNESNEWTONTR`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESVIGetInactiveSet()`, `DMSetVI()`, `SNESVISetRedundancyCheck()`
751: M*/
752: PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONRSLS(SNES snes)
753: {
754:   SNES_VINEWTONRSLS *vi;
755:   SNESLineSearch     linesearch;

757:   PetscFunctionBegin;
758:   snes->ops->reset          = SNESReset_VINEWTONRSLS;
759:   snes->ops->setup          = SNESSetUp_VINEWTONRSLS;
760:   snes->ops->solve          = SNESSolve_VINEWTONRSLS;
761:   snes->ops->destroy        = SNESDestroy_VI;
762:   snes->ops->setfromoptions = SNESSetFromOptions_VI;
763:   snes->ops->view           = NULL;
764:   snes->ops->converged      = SNESConvergedDefault_VI;

766:   snes->usesksp = PETSC_TRUE;
767:   snes->usesnpc = PETSC_FALSE;

769:   PetscCall(SNESGetLineSearch(snes, &linesearch));
770:   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
771:   PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0));

773:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

775:   PetscCall(SNESParametersInitialize(snes));

777:   PetscCall(PetscNew(&vi));
778:   snes->data          = (void *)vi;
779:   vi->checkredundancy = NULL;

781:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI));
782:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI));
783:   PetscFunctionReturn(PETSC_SUCCESS);
784: }