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 parameter:
  9: +  snes - the `SNES` context
 10: -  compute - function that computes the bounds

 12:  Calling Sequence of function:
 13:   PetscErrorCode compute(SNES snes,Vec lower,Vec higher, void *ctx)

 15: + snes - the `SNES` context
 16: . lower - vector to hold lower bounds
 17: - higher - vector to hold upper bounds

 19:    Level: advanced

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

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

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

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

 31: .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `DMSetVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`,
 32:           'SNESSetType()`
 33: @*/
 34: PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES, Vec, Vec))
 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, SNESVIComputeVariableBoundsFunction compute)
 47: {
 48:   PetscFunctionBegin;
 49:   snes->ops->computevariablebounds = compute;
 50:   PetscFunctionReturn(PETSC_SUCCESS);
 51: }

 53: PetscErrorCode SNESVIMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy)
 54: {
 55:   Vec         X, F, Finactive;
 56:   IS          isactive;
 57:   PetscViewer viewer = (PetscViewer)dummy;

 59:   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(VecCopy(F, Finactive));
 65:   PetscCall(VecISSet(Finactive, isactive, 0.0));
 66:   PetscCall(ISDestroy(&isactive));
 67:   PetscCall(VecView(Finactive, viewer));
 68:   PetscCall(VecDestroy(&Finactive));
 69:   PetscFunctionReturn(PETSC_SUCCESS);
 70: }

 72: PetscErrorCode SNESMonitorVI(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy)
 73: {
 74:   PetscViewer        viewer = (PetscViewer)dummy;
 75:   const PetscScalar *x, *xl, *xu, *f;
 76:   PetscInt           i, n, act[2] = {0, 0}, fact[2], N;
 77:   /* Number of components that actually hit the bounds (c.f. active variables) */
 78:   PetscInt  act_bound[2] = {0, 0}, fact_bound[2];
 79:   PetscReal rnorm, fnorm, zerotolerance = snes->vizerotolerance;
 80:   double    tmp;

 82:   PetscFunctionBegin;
 84:   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
 85:   PetscCall(VecGetSize(snes->vec_sol, &N));
 86:   PetscCall(VecGetArrayRead(snes->xl, &xl));
 87:   PetscCall(VecGetArrayRead(snes->xu, &xu));
 88:   PetscCall(VecGetArrayRead(snes->vec_sol, &x));
 89:   PetscCall(VecGetArrayRead(snes->vec_func, &f));

 91:   rnorm = 0.0;
 92:   for (i = 0; i < n; i++) {
 93:     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]);
 94:     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++;
 95:     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++;
 96:     else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Can never get here");
 97:   }

 99:   for (i = 0; i < n; i++) {
100:     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++;
101:     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++;
102:   }
103:   PetscCall(VecRestoreArrayRead(snes->vec_func, &f));
104:   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
105:   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
106:   PetscCall(VecRestoreArrayRead(snes->vec_sol, &x));
107:   PetscCall(MPIU_Allreduce(&rnorm, &fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
108:   PetscCall(MPIU_Allreduce(act, fact, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
109:   PetscCall(MPIU_Allreduce(act_bound, fact_bound, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
110:   fnorm = PetscSqrtReal(fnorm);

112:   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
113:   if (snes->ntruebounds) tmp = ((double)(fact[0] + fact[1])) / ((double)snes->ntruebounds);
114:   else tmp = 0.0;
115:   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));

117:   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
118:   PetscFunctionReturn(PETSC_SUCCESS);
119: }

121: /*
122:      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
123:     || 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
124:     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
125:     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
126: */
127: PetscErrorCode SNESVICheckLocalMin_Private(SNES snes, Mat A, Vec F, Vec W, PetscReal fnorm, PetscBool *ismin)
128: {
129:   PetscReal a1;
130:   PetscBool hastranspose;

132:   PetscFunctionBegin;
133:   *ismin = PETSC_FALSE;
134:   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
135:   if (hastranspose) {
136:     /* Compute || J^T F|| */
137:     PetscCall(MatMultTranspose(A, F, W));
138:     PetscCall(VecNorm(W, NORM_2, &a1));
139:     PetscCall(PetscInfo(snes, "|| J^T F|| %g near zero implies found a local minimum\n", (double)(a1 / fnorm)));
140:     if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE;
141:   } else {
142:     Vec         work;
143:     PetscScalar result;
144:     PetscReal   wnorm;

146:     PetscCall(VecSetRandom(W, NULL));
147:     PetscCall(VecNorm(W, NORM_2, &wnorm));
148:     PetscCall(VecDuplicate(W, &work));
149:     PetscCall(MatMult(A, W, work));
150:     PetscCall(VecDot(F, work, &result));
151:     PetscCall(VecDestroy(&work));
152:     a1 = PetscAbsScalar(result) / (fnorm * wnorm);
153:     PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n", (double)a1));
154:     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
155:   }
156:   PetscFunctionReturn(PETSC_SUCCESS);
157: }

159: /*
160:      Checks if J^T(F - J*X) = 0
161: */
162: PetscErrorCode SNESVICheckResidual_Private(SNES snes, Mat A, Vec F, Vec X, Vec W1, Vec W2)
163: {
164:   PetscReal a1, a2;
165:   PetscBool hastranspose;

167:   PetscFunctionBegin;
168:   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
169:   if (hastranspose) {
170:     PetscCall(MatMult(A, X, W1));
171:     PetscCall(VecAXPY(W1, -1.0, F));

173:     /* Compute || J^T W|| */
174:     PetscCall(MatMultTranspose(A, W1, W2));
175:     PetscCall(VecNorm(W1, NORM_2, &a1));
176:     PetscCall(VecNorm(W2, NORM_2, &a2));
177:     if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n", (double)(a2 / a1)));
178:   }
179:   PetscFunctionReturn(PETSC_SUCCESS);
180: }

182: /*
183:   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.

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

196:   *reason = SNES_CONVERGED_ITERATING;

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

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

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

225:    Input Parameters:
226: .  SNES - nonlinear solver context

228:    Output Parameters:
229: .  X - Bound projected X

231: */

233: PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X)
234: {
235:   const PetscScalar *xl, *xu;
236:   PetscScalar       *x;
237:   PetscInt           i, n;

239:   PetscFunctionBegin;
240:   PetscCall(VecGetLocalSize(X, &n));
241:   PetscCall(VecGetArray(X, &x));
242:   PetscCall(VecGetArrayRead(snes->xl, &xl));
243:   PetscCall(VecGetArrayRead(snes->xu, &xu));

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

255: /*
256:    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables

258:    Input parameter:
259: .  snes - the SNES context
260: .  X    - the snes solution vector
261: .  F    - the nonlinear function vector

263:    Output parameter:
264: .  ISact - active set index set
265:  */
266: PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact)
267: {
268:   Vec                Xl = snes->xl, Xu = snes->xu;
269:   const PetscScalar *x, *f, *xl, *xu;
270:   PetscInt          *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0;
271:   PetscReal          zerotolerance = snes->vizerotolerance;

273:   PetscFunctionBegin;
274:   PetscCall(VecGetLocalSize(X, &nlocal));
275:   PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh));
276:   PetscCall(VecGetArrayRead(X, &x));
277:   PetscCall(VecGetArrayRead(Xl, &xl));
278:   PetscCall(VecGetArrayRead(Xu, &xu));
279:   PetscCall(VecGetArrayRead(F, &f));
280:   /* Compute active set size */
281:   for (i = 0; i < nlocal; i++) {
282:     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++;
283:   }

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

287:   /* Set active set indices */
288:   for (i = 0; i < nlocal; i++) {
289:     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;
290:   }

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

295:   PetscCall(VecRestoreArrayRead(X, &x));
296:   PetscCall(VecRestoreArrayRead(Xl, &xl));
297:   PetscCall(VecRestoreArrayRead(Xu, &xu));
298:   PetscCall(VecRestoreArrayRead(F, &f));
299:   PetscFunctionReturn(PETSC_SUCCESS);
300: }

302: PetscErrorCode SNESVICreateIndexSets_RS(SNES snes, Vec X, Vec F, IS *ISact, IS *ISinact)
303: {
304:   PetscInt rstart, rend;

306:   PetscFunctionBegin;
307:   PetscCall(SNESVIGetActiveSetIS(snes, X, F, ISact));
308:   PetscCall(VecGetOwnershipRange(X, &rstart, &rend));
309:   PetscCall(ISComplement(*ISact, rstart, rend, ISinact));
310:   PetscFunctionReturn(PETSC_SUCCESS);
311: }

313: PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm)
314: {
315:   const PetscScalar *x, *xl, *xu, *f;
316:   PetscInt           i, n;
317:   PetscReal          rnorm, zerotolerance = snes->vizerotolerance;

319:   PetscFunctionBegin;
320:   PetscCall(VecGetLocalSize(X, &n));
321:   PetscCall(VecGetArrayRead(snes->xl, &xl));
322:   PetscCall(VecGetArrayRead(snes->xu, &xu));
323:   PetscCall(VecGetArrayRead(X, &x));
324:   PetscCall(VecGetArrayRead(F, &f));
325:   rnorm = 0.0;
326:   for (i = 0; i < n; i++) {
327:     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]);
328:   }
329:   PetscCall(VecRestoreArrayRead(F, &f));
330:   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
331:   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
332:   PetscCall(VecRestoreArrayRead(X, &x));
333:   PetscCall(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
334:   *fnorm = PetscSqrtReal(*fnorm);
335:   PetscFunctionReturn(PETSC_SUCCESS);
336: }

338: PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu)
339: {
340:   PetscFunctionBegin;
341:   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
342:   PetscFunctionReturn(PETSC_SUCCESS);
343: }

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

349:    Input Parameter:
350: .  snes - the SNES context

352:    Application Interface Routine: SNESSetUp()

354:    Notes:
355:    For basic use of the SNES solvers, the user need not explicitly call
356:    SNESSetUp(), since these actions will automatically occur during
357:    the call to SNESSolve().
358:  */
359: PetscErrorCode SNESSetUp_VI(SNES snes)
360: {
361:   PetscInt i_start[3], i_end[3];

363:   PetscFunctionBegin;
364:   PetscCall(SNESSetWorkVecs(snes, 1));
365:   PetscCall(SNESSetUpMatrices(snes));

367:   if (!snes->ops->computevariablebounds && snes->dm) {
368:     PetscBool flag;
369:     PetscCall(DMHasVariableBounds(snes->dm, &flag));
370:     if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
371:   }
372:   if (!snes->usersetbounds) {
373:     if (snes->ops->computevariablebounds) {
374:       if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
375:       if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
376:       PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
377:     } else if (!snes->xl && !snes->xu) {
378:       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
379:       PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
380:       PetscCall(VecSet(snes->xl, PETSC_NINFINITY));
381:       PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
382:       PetscCall(VecSet(snes->xu, PETSC_INFINITY));
383:     } else {
384:       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
385:       PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end));
386:       PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1));
387:       PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2));
388:       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]))
389:         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.");
390:     }
391:   }
392:   PetscFunctionReturn(PETSC_SUCCESS);
393: }
394: PetscErrorCode SNESReset_VI(SNES snes)
395: {
396:   PetscFunctionBegin;
397:   PetscCall(VecDestroy(&snes->xl));
398:   PetscCall(VecDestroy(&snes->xu));
399:   snes->usersetbounds = PETSC_FALSE;
400:   PetscFunctionReturn(PETSC_SUCCESS);
401: }

403: /*
404:    SNESDestroy_VI - Destroys the private SNES_VI context that was created
405:    with SNESCreate_VI().

407:    Input Parameter:
408: .  snes - the SNES context

410:    Application Interface Routine: SNESDestroy()
411:  */
412: PetscErrorCode SNESDestroy_VI(SNES snes)
413: {
414:   PetscFunctionBegin;
415:   PetscCall(PetscFree(snes->data));

417:   /* clear composed functions */
418:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL));
419:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL));
420:   PetscFunctionReturn(PETSC_SUCCESS);
421: }

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

427:    Input Parameters:
428: +  snes - the `SNES` context.
429: .  xl   - lower bound.
430: -  xu   - upper bound.

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

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

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

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

443:    Level: advanced

445: .seealso: [](sec_vi), `SNES`, `SNESVIGetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, 'SNESSetType()`
446: @*/
447: PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
448: {
449:   PetscErrorCode (*f)(SNES, Vec, Vec);

451:   PetscFunctionBegin;
455:   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f));
456:   if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu));
457:   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
458:   snes->usersetbounds = PETSC_TRUE;
459:   PetscFunctionReturn(PETSC_SUCCESS);
460: }

462: PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu)
463: {
464:   const PetscScalar *xxl, *xxu;
465:   PetscInt           i, n, cnt = 0;

467:   PetscFunctionBegin;
468:   PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL));
469:   PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first");
470:   {
471:     PetscInt xlN, xuN, N;
472:     PetscCall(VecGetSize(xl, &xlN));
473:     PetscCall(VecGetSize(xu, &xuN));
474:     PetscCall(VecGetSize(snes->vec_func, &N));
475:     PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N);
476:     PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N);
477:   }
478:   PetscCall(PetscObjectReference((PetscObject)xl));
479:   PetscCall(PetscObjectReference((PetscObject)xu));
480:   PetscCall(VecDestroy(&snes->xl));
481:   PetscCall(VecDestroy(&snes->xu));
482:   snes->xl = xl;
483:   snes->xu = xu;
484:   PetscCall(VecGetLocalSize(xl, &n));
485:   PetscCall(VecGetArrayRead(xl, &xxl));
486:   PetscCall(VecGetArrayRead(xu, &xxu));
487:   for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));

489:   PetscCall(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
490:   PetscCall(VecRestoreArrayRead(xl, &xxl));
491:   PetscCall(VecRestoreArrayRead(xu, &xxu));
492:   PetscFunctionReturn(PETSC_SUCCESS);
493: }

495: /*@
496:    SNESVIGetVariableBounds - Gets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving
497:    (differential) variable inequalities.

499:    Input Parameters:
500: +  snes - the `SNES` context.
501: .  xl   - lower bound (may be `NULL`)
502: -  xu   - upper bound (may be `NULL`)

504:    Notes:
505:    These vectors are owned by the `SNESVI` and should not be destroyed by the caller

507:    Level: advanced

509: .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, SNESVINEWTONRSLS, SNESVINEWTONSSLS, 'SNESSetType()`
510: @*/
511: PetscErrorCode SNESVIGetVariableBounds(SNES snes, Vec *xl, Vec *xu)
512: {
513:   PetscFunctionBegin;
514:   PetscCheck(snes->usersetbounds, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must set SNESVI bounds before calling SNESVIGetVariableBounds()");
515:   if (xl) *xl = snes->xl;
516:   if (xu) *xu = snes->xu;
517:   PetscFunctionReturn(PETSC_SUCCESS);
518: }

520: PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject)
521: {
522:   PetscBool flg = PETSC_FALSE;

524:   PetscFunctionBegin;
525:   PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options");
526:   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL));
527:   PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL));
528:   if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL));
529:   flg = PETSC_FALSE;
530:   PetscCall(PetscOptionsBool("-snes_vi_monitor_residual", "Monitor residual all non-active variables; using zero for active constraints", "SNESMonitorVIResidual", flg, &flg, NULL));
531:   if (flg) PetscCall(SNESMonitorSet(snes, SNESVIMonitorResidual, PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)), NULL));
532:   PetscOptionsHeadEnd();
533:   PetscFunctionReturn(PETSC_SUCCESS);
534: }