Actual source code: vi.c

  1: #include <petsc/private/snesimpl.h>
  2: #include <petscdm.h>

  4: /*@C
  5:   SNESVISetComputeVariableBounds - Sets a function that is called to compute the bounds on variable for
  6:   (differential) variable inequalities.

  8:   Input Parameters:
  9: + snes    - the `SNES` context
 10: - compute - function that computes the bounds

 12:   Calling sequence of `compute`:
 13: + snes   - the `SNES` context
 14: . lower  - vector to hold lower bounds
 15: - higher - vector to hold upper bounds

 17:   Level: advanced

 19:   Notes:
 20:   Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.

 22:   For entries with no bounds you can set `PETSC_NINFINITY` or `PETSC_INFINITY`

 24:   You may use `SNESVISetVariableBounds()` to provide the bounds once if they will never change

 26:   If you have associated a `DM` with the `SNES` and provided a function to the `DM` via `DMSetVariableBounds()` that will be used automatically
 27:   to provide the bounds and you need not use this function.

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

 31: .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `DMSetVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`,
 32:           `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
 33: @*/
 34: PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES snes, Vec lower, Vec higher))
 35: {
 36:   PetscErrorCode (*f)(SNES, PetscErrorCode (*)(SNES, Vec, Vec));

 38:   PetscFunctionBegin;
 40:   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", &f));
 41:   if (f) PetscUseMethod(snes, "SNESVISetComputeVariableBounds_C", (SNES, PetscErrorCode (*)(SNES, Vec, Vec)), (snes, compute));
 42:   else PetscCall(SNESVISetComputeVariableBounds_VI(snes, compute));
 43:   PetscFunctionReturn(PETSC_SUCCESS);
 44: }

 46: PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes, SNESVIComputeVariableBoundsFn *compute)
 47: {
 48:   PetscFunctionBegin;
 49:   snes->ops->computevariablebounds = compute;
 50:   PetscFunctionReturn(PETSC_SUCCESS);
 51: }

 53: static PetscErrorCode SNESVIMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
 54: {
 55:   Vec X, F, Finactive;
 56:   IS  isactive;

 58:   PetscFunctionBegin;
 60:   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
 61:   PetscCall(SNESGetSolution(snes, &X));
 62:   PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive));
 63:   PetscCall(VecDuplicate(F, &Finactive));
 64:   PetscCall(PetscObjectCompose((PetscObject)Finactive, "__Vec_bc_zero__", (PetscObject)snes));
 65:   PetscCall(VecCopy(F, Finactive));
 66:   PetscCall(VecISSet(Finactive, isactive, 0.0));
 67:   PetscCall(ISDestroy(&isactive));
 68:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
 69:   PetscCall(VecView(Finactive, vf->viewer));
 70:   PetscCall(PetscViewerPopFormat(vf->viewer));
 71:   PetscCall(PetscObjectCompose((PetscObject)Finactive, "__Vec_bc_zero__", NULL));
 72:   PetscCall(VecDestroy(&Finactive));
 73:   PetscFunctionReturn(PETSC_SUCCESS);
 74: }

 76: static PetscErrorCode SNESVIMonitorActive(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
 77: {
 78:   Vec X, F, A;
 79:   IS  isactive;

 81:   PetscFunctionBegin;
 83:   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
 84:   PetscCall(SNESGetSolution(snes, &X));
 85:   PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive));
 86:   PetscCall(VecDuplicate(F, &A));
 87:   PetscCall(PetscObjectCompose((PetscObject)A, "__Vec_bc_zero__", (PetscObject)snes));
 88:   PetscCall(VecISSet(A, isactive, 1.));
 89:   PetscCall(ISDestroy(&isactive));
 90:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
 91:   PetscCall(VecView(A, vf->viewer));
 92:   PetscCall(PetscViewerPopFormat(vf->viewer));
 93:   PetscCall(PetscObjectCompose((PetscObject)A, "__Vec_bc_zero__", NULL));
 94:   PetscCall(VecDestroy(&A));
 95:   PetscFunctionReturn(PETSC_SUCCESS);
 96: }

 98: static PetscErrorCode SNESMonitorVI(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy)
 99: {
100:   PetscViewer        viewer = (PetscViewer)dummy;
101:   const PetscScalar *x, *xl, *xu, *f;
102:   PetscInt           i, n, act[2] = {0, 0}, fact[2], N;
103:   /* Number of components that actually hit the bounds (c.f. active variables) */
104:   PetscInt  act_bound[2] = {0, 0}, fact_bound[2];
105:   PetscReal rnorm, fnorm, zerotolerance = snes->vizerotolerance;
106:   double    tmp;

108:   PetscFunctionBegin;
110:   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
111:   PetscCall(VecGetSize(snes->vec_sol, &N));
112:   PetscCall(VecGetArrayRead(snes->xl, &xl));
113:   PetscCall(VecGetArrayRead(snes->xu, &xu));
114:   PetscCall(VecGetArrayRead(snes->vec_sol, &x));
115:   PetscCall(VecGetArrayRead(snes->vec_func, &f));

117:   rnorm = 0.0;
118:   for (i = 0; i < n; i++) {
119:     if ((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0)) rnorm += PetscRealPart(PetscConj(f[i]) * f[i]);
120:     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++;
121:     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++;
122:     else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Can never get here");
123:   }

125:   for (i = 0; i < n; i++) {
126:     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++;
127:     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++;
128:   }
129:   PetscCall(VecRestoreArrayRead(snes->vec_func, &f));
130:   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
131:   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
132:   PetscCall(VecRestoreArrayRead(snes->vec_sol, &x));
133:   PetscCallMPI(MPIU_Allreduce(&rnorm, &fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
134:   PetscCallMPI(MPIU_Allreduce(act, fact, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
135:   PetscCallMPI(MPIU_Allreduce(act_bound, fact_bound, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
136:   fnorm = PetscSqrtReal(fnorm);

138:   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
139:   if (snes->ntruebounds) tmp = ((double)(fact[0] + fact[1])) / ((double)snes->ntruebounds);
140:   else tmp = 0.0;
141:   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES VI Function norm %g Active lower constraints %" PetscInt_FMT "/%" PetscInt_FMT " upper constraints %" PetscInt_FMT "/%" PetscInt_FMT " Percent of total %g Percent of bounded %g\n", its, (double)fnorm, fact[0], fact_bound[0], fact[1], fact_bound[1], ((double)(fact[0] + fact[1])) / ((double)N), tmp));

143:   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
144:   PetscFunctionReturn(PETSC_SUCCESS);
145: }

147: /*
148:      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
149:     || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that
150:     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
151:     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
152: */
153: PetscErrorCode SNESVICheckLocalMin_Private(SNES snes, Mat A, Vec F, Vec W, PetscReal fnorm, PetscBool *ismin)
154: {
155:   PetscReal a1;
156:   PetscBool hastranspose;

158:   PetscFunctionBegin;
159:   *ismin = PETSC_FALSE;
160:   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
161:   if (hastranspose) {
162:     /* Compute || J^T F|| */
163:     PetscCall(MatMultTranspose(A, F, W));
164:     PetscCall(VecNorm(W, NORM_2, &a1));
165:     PetscCall(PetscInfo(snes, "|| J^T F|| %g near zero implies found a local minimum\n", (double)(a1 / fnorm)));
166:     if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE;
167:   } else {
168:     Vec         work;
169:     PetscScalar result;
170:     PetscReal   wnorm;

172:     PetscCall(VecSetRandom(W, NULL));
173:     PetscCall(VecNorm(W, NORM_2, &wnorm));
174:     PetscCall(VecDuplicate(W, &work));
175:     PetscCall(MatMult(A, W, work));
176:     PetscCall(VecDot(F, work, &result));
177:     PetscCall(VecDestroy(&work));
178:     a1 = PetscAbsScalar(result) / (fnorm * wnorm);
179:     PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n", (double)a1));
180:     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
181:   }
182:   PetscFunctionReturn(PETSC_SUCCESS);
183: }

185: /*
186:   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.

188:   Notes:
189:   The convergence criterion currently implemented is
190:   merit < abstol
191:   merit < rtol*merit_initial
192: */
193: PetscErrorCode SNESConvergedDefault_VI(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gradnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
194: {
195:   PetscFunctionBegin;
197:   PetscAssertPointer(reason, 6);

199:   *reason = SNES_CONVERGED_ITERATING;

201:   if (!it) {
202:     /* set parameter for default relative tolerance convergence test */
203:     snes->ttol = fnorm * snes->rtol;
204:   }
205:   if (fnorm != fnorm) {
206:     PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n"));
207:     *reason = SNES_DIVERGED_FUNCTION_NANORINF;
208:   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
209:     PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g\n", (double)fnorm, (double)snes->abstol));
210:     *reason = SNES_CONVERGED_FNORM_ABS;
211:   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
212:     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs));
213:     *reason = SNES_DIVERGED_FUNCTION_COUNT;
214:   }

216:   if (it && !*reason) {
217:     if (fnorm < snes->ttol) {
218:       PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g (relative tolerance)\n", (double)fnorm, (double)snes->ttol));
219:       *reason = SNES_CONVERGED_FNORM_RELATIVE;
220:     }
221:   }
222:   PetscFunctionReturn(PETSC_SUCCESS);
223: }

225: /*
226:    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.

228:    Input Parameters:
229: .  SNES - nonlinear solver context

231:    Output Parameters:
232: .  X - Bound projected X

234: */

236: PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X)
237: {
238:   const PetscScalar *xl, *xu;
239:   PetscScalar       *x;
240:   PetscInt           i, n;

242:   PetscFunctionBegin;
243:   PetscCall(VecGetLocalSize(X, &n));
244:   PetscCall(VecGetArray(X, &x));
245:   PetscCall(VecGetArrayRead(snes->xl, &xl));
246:   PetscCall(VecGetArrayRead(snes->xu, &xu));

248:   for (i = 0; i < n; i++) {
249:     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
250:     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
251:   }
252:   PetscCall(VecRestoreArray(X, &x));
253:   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
254:   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
255:   PetscFunctionReturn(PETSC_SUCCESS);
256: }

258: /*@
259:   SNESVIGetActiveSetIS - Gets the global indices for the active set variables

261:   Input Parameters:
262: + snes - the `SNES` context
263: . X    - the `snes` solution vector
264: - F    - the nonlinear function vector

266:   Output Parameter:
267: . ISact - active set index set

269:   Level: developer

271:   Note:
272:   See `SNESVINEWTONRSLS` for a concise description of the active and inactive sets

274: .seealso: [](ch_snes), `SNES`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`
275: @*/
276: PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact)
277: {
278:   Vec                Xl = snes->xl, Xu = snes->xu;
279:   const PetscScalar *x, *f, *xl, *xu;
280:   PetscInt          *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0;
281:   PetscReal          zerotolerance = snes->vizerotolerance;

283:   PetscFunctionBegin;
284:   PetscCall(VecGetLocalSize(X, &nlocal));
285:   PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh));
286:   PetscCall(VecGetArrayRead(X, &x));
287:   PetscCall(VecGetArrayRead(Xl, &xl));
288:   PetscCall(VecGetArrayRead(Xu, &xu));
289:   PetscCall(VecGetArrayRead(F, &f));
290:   /* Compute active set size */
291:   for (i = 0; i < nlocal; i++) {
292:     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) nloc_isact++;
293:   }

295:   PetscCall(PetscMalloc1(nloc_isact, &idx_act));

297:   /* Set active set indices */
298:   for (i = 0; i < nlocal; i++) {
299:     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) idx_act[i1++] = ilow + i;
300:   }

302:   /* Create active set IS */
303:   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact));

305:   PetscCall(VecRestoreArrayRead(X, &x));
306:   PetscCall(VecRestoreArrayRead(Xl, &xl));
307:   PetscCall(VecRestoreArrayRead(Xu, &xu));
308:   PetscCall(VecRestoreArrayRead(F, &f));
309:   PetscFunctionReturn(PETSC_SUCCESS);
310: }

312: /*@
313:   SNESVIComputeInactiveSetFnorm - Computes the function norm for variational inequalities on the inactive set

315:   Input Parameters:
316: + snes - the `SNES` context
317: . F    - the nonlinear function vector
318: - X    - the `SNES` solution vector

320:   Output Parameter:
321: . fnorm - the function norm

323:   Level: developer

325:   Note:
326:   See `SNESVINEWTONRSLS` for a concise description of the active and inactive sets

328: .seealso: [](ch_snes), `SNES`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `SNESLineSearchSetVIFunctions()`
329: @*/
330: PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm)
331: {
332:   const PetscScalar *x, *xl, *xu, *f;
333:   PetscInt           i, n;
334:   PetscReal          zerotolerance = snes->vizerotolerance;

336:   PetscFunctionBegin;
338:   PetscAssertPointer(fnorm, 4);
339:   PetscCall(VecGetLocalSize(X, &n));
340:   PetscCall(VecGetArrayRead(snes->xl, &xl));
341:   PetscCall(VecGetArrayRead(snes->xu, &xu));
342:   PetscCall(VecGetArrayRead(X, &x));
343:   PetscCall(VecGetArrayRead(F, &f));
344:   *fnorm = 0.0;
345:   for (i = 0; i < n; i++) {
346:     if ((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0)) *fnorm += PetscRealPart(PetscConj(f[i]) * f[i]);
347:   }
348:   PetscCall(VecRestoreArrayRead(F, &f));
349:   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
350:   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
351:   PetscCall(VecRestoreArrayRead(X, &x));
352:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
353:   *fnorm = PetscSqrtReal(*fnorm);
354:   PetscFunctionReturn(PETSC_SUCCESS);
355: }

357: /*@
358:   SNESVIComputeInactiveSetFtY - Computes the directional derivative for variational inequalities on the inactive set,
359:   assuming that there exists some $G(x)$ for which the `SNESFunctionFn` $F(x) = grad G(x)$ (relevant for some line search algorithms)

361:   Input Parameters:
362: + snes - the `SNES` context
363: . F    - the nonlinear function vector
364: . X    - the `SNES` solution vector
365: - Y    - the direction vector

367:   Output Parameter:
368: . fty - the directional derivative

370:   Level: developer

372:   Note:
373:   See `SNESVINEWTONRSLS` for a concise description of the active and inactive sets

375: .seealso: [](ch_snes), `SNES`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`
376: @*/
377: PetscErrorCode SNESVIComputeInactiveSetFtY(SNES snes, Vec F, Vec X, Vec Y, PetscScalar *fty)
378: {
379:   const PetscScalar *x, *xl, *xu, *y, *f;
380:   PetscInt           i, n;
381:   PetscReal          zerotolerance = snes->vizerotolerance;

383:   PetscFunctionBegin;
385:   PetscAssertPointer(fty, 5);
386:   PetscCall(VecGetLocalSize(X, &n));
387:   PetscCall(VecGetArrayRead(F, &f));
388:   PetscCall(VecGetArrayRead(X, &x));
389:   PetscCall(VecGetArrayRead(snes->xl, &xl));
390:   PetscCall(VecGetArrayRead(snes->xu, &xu));
391:   PetscCall(VecGetArrayRead(Y, &y));
392:   *fty = 0.0;
393:   for (i = 0; i < n; i++) {
394:     if ((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0)) *fty += f[i] * PetscConj(y[i]);
395:   }
396:   PetscCall(VecRestoreArrayRead(F, &f));
397:   PetscCall(VecRestoreArrayRead(X, &x));
398:   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
399:   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
400:   PetscCall(VecRestoreArrayRead(Y, &y));
401:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, fty, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
402:   PetscFunctionReturn(PETSC_SUCCESS);
403: }

405: static PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu)
406: {
407:   PetscFunctionBegin;
408:   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
409:   PetscFunctionReturn(PETSC_SUCCESS);
410: }

412: /*
413:    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
414:    of the SNESVI nonlinear solver.

416:    Input Parameter:
417: .  snes - the SNES context

419:    Application Interface Routine: SNESSetUp()

421:    Notes:
422:    For basic use of the SNES solvers, the user need not explicitly call
423:    SNESSetUp(), since these actions will automatically occur during
424:    the call to SNESSolve().
425:  */
426: PetscErrorCode SNESSetUp_VI(SNES snes)
427: {
428:   PetscInt i_start[3], i_end[3];

430:   PetscFunctionBegin;
431:   PetscCall(SNESSetWorkVecs(snes, 1));
432:   PetscCall(SNESSetUpMatrices(snes));

434:   if (!snes->ops->computevariablebounds && snes->dm) {
435:     PetscBool flag;
436:     PetscCall(DMHasVariableBounds(snes->dm, &flag));
437:     if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
438:   }
439:   if (!snes->usersetbounds) {
440:     if (snes->ops->computevariablebounds) {
441:       if (!snes->xl) PetscCall(VecDuplicate(snes->work[0], &snes->xl));
442:       if (!snes->xu) PetscCall(VecDuplicate(snes->work[0], &snes->xu));
443:       PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
444:     } else if (!snes->xl && !snes->xu) {
445:       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
446:       PetscCall(VecDuplicate(snes->work[0], &snes->xl));
447:       PetscCall(VecSet(snes->xl, PETSC_NINFINITY));
448:       PetscCall(VecDuplicate(snes->work[0], &snes->xu));
449:       PetscCall(VecSet(snes->xu, PETSC_INFINITY));
450:     } else {
451:       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
452:       PetscCall(VecGetOwnershipRange(snes->work[0], i_start, i_end));
453:       PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1));
454:       PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2));
455:       if ((i_start[0] != i_start[1]) || (i_start[0] != i_start[2]) || (i_end[0] != i_end[1]) || (i_end[0] != i_end[2]))
456:         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Distribution of lower bound, upper bound and the solution vector should be identical across all the processors.");
457:     }
458:   }
459:   PetscFunctionReturn(PETSC_SUCCESS);
460: }
461: PetscErrorCode SNESReset_VI(SNES snes)
462: {
463:   PetscFunctionBegin;
464:   PetscCall(VecDestroy(&snes->xl));
465:   PetscCall(VecDestroy(&snes->xu));
466:   snes->usersetbounds = PETSC_FALSE;
467:   PetscFunctionReturn(PETSC_SUCCESS);
468: }

470: /*
471:    SNESDestroy_VI - Destroys the private SNES_VI context that was created
472:    with SNESCreate_VI().

474:    Input Parameter:
475: .  snes - the SNES context

477:    Application Interface Routine: SNESDestroy()
478:  */
479: PetscErrorCode SNESDestroy_VI(SNES snes)
480: {
481:   PetscFunctionBegin;
482:   PetscCall(PetscFree(snes->data));

484:   /* clear composed functions */
485:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL));
486:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL));
487:   PetscFunctionReturn(PETSC_SUCCESS);
488: }

490: /*@
491:   SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. `xl` <= x <= `xu`. This allows solving
492:   (differential) variable inequalities.

494:   Input Parameters:
495: + snes - the `SNES` context.
496: . xl   - lower bound.
497: - xu   - upper bound.

499:   Level: advanced

501:   Notes:
502:   If this routine is not called then the lower and upper bounds are set to
503:   `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`.

505:   Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS` or semi-smooth `SNESVINEWTONSSLS` solvers.

507:   For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY`

509:   `SNESVISetComputeVariableBounds()` can be used to provide a function that computes the bounds. This should be used if you are using, for example, grid
510:   sequencing and need bounds set for a variety of vectors

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

514: .seealso: [](sec_vi), `SNES`, `SNESVIGetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
515: @*/
516: PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
517: {
518:   PetscErrorCode (*f)(SNES, Vec, Vec);

520:   PetscFunctionBegin;
524:   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f));
525:   if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu));
526:   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
527:   snes->usersetbounds = PETSC_TRUE;
528:   PetscFunctionReturn(PETSC_SUCCESS);
529: }

531: PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu)
532: {
533:   const PetscScalar *xxl, *xxu;
534:   PetscInt           i, n, cnt = 0;

536:   PetscFunctionBegin;
537:   PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL));
538:   PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first");
539:   {
540:     PetscInt xlN, xuN, N;
541:     PetscCall(VecGetSize(xl, &xlN));
542:     PetscCall(VecGetSize(xu, &xuN));
543:     PetscCall(VecGetSize(snes->vec_func, &N));
544:     PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N);
545:     PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N);
546:   }
547:   PetscCall(PetscObjectReference((PetscObject)xl));
548:   PetscCall(PetscObjectReference((PetscObject)xu));
549:   PetscCall(VecDestroy(&snes->xl));
550:   PetscCall(VecDestroy(&snes->xu));
551:   snes->xl = xl;
552:   snes->xu = xu;
553:   PetscCall(VecGetLocalSize(xl, &n));
554:   PetscCall(VecGetArrayRead(xl, &xxl));
555:   PetscCall(VecGetArrayRead(xu, &xxu));
556:   for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));

558:   PetscCallMPI(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
559:   PetscCall(VecRestoreArrayRead(xl, &xxl));
560:   PetscCall(VecRestoreArrayRead(xu, &xxu));
561:   PetscFunctionReturn(PETSC_SUCCESS);
562: }

564: /*@
565:   SNESVIGetVariableBounds - Gets the lower and upper bounds for the solution vector. `xl` <= x <= `xu`. These are used in solving
566:   (differential) variable inequalities.

568:   Input Parameters:
569: + snes - the `SNES` context.
570: . xl   - lower bound (may be `NULL`)
571: - xu   - upper bound (may be `NULL`)

573:   Level: advanced

575:   Note:
576:   These vectors are owned by the `SNESVI` and should not be destroyed by the caller

578: .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
579: @*/
580: PetscErrorCode SNESVIGetVariableBounds(SNES snes, Vec *xl, Vec *xu)
581: {
582:   PetscFunctionBegin;
583:   PetscCheck(snes->usersetbounds, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must set SNESVI bounds before calling SNESVIGetVariableBounds()");
584:   if (xl) *xl = snes->xl;
585:   if (xu) *xu = snes->xu;
586:   PetscFunctionReturn(PETSC_SUCCESS);
587: }

589: PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems PetscOptionsObject)
590: {
591:   PetscBool flg = PETSC_FALSE;

593:   PetscFunctionBegin;
594:   PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options");
595:   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL));
596:   PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL));
597:   if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL));
598:   flg = PETSC_FALSE;
599:   PetscCall(SNESMonitorSetFromOptions(snes, "-snes_vi_monitor_residual", "View residual at each iteration, using zero for active constraints", "SNESVIMonitorResidual", SNESVIMonitorResidual, NULL));
600:   PetscCall(SNESMonitorSetFromOptions(snes, "-snes_vi_monitor_active", "View active set at each iteration, using zero for inactive dofs", "SNESVIMonitorActive", SNESVIMonitorActive, NULL));
601:   PetscOptionsHeadEnd();
602:   PetscFunctionReturn(PETSC_SUCCESS);
603: }