xref: /petsc/src/snes/impls/vi/vi.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2 #include <petscdm.h>
3 
4 /*@C
5   SNESVISetComputeVariableBounds - Sets a function that is called to compute the bounds on variable for
6   (differential) variable inequalities.
7 
8   Input Parameters:
9 + snes    - the `SNES` context
10 - compute - function that computes the bounds
11 
12   Calling sequence of `compute`:
13 + snes   - the `SNES` context
14 . lower  - vector to hold lower bounds
15 - higher - vector to hold upper bounds
16 
17   Level: advanced
18 
19   Notes:
20   Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.
21 
22   For entries with no bounds you can set `PETSC_NINFINITY` or `PETSC_INFINITY`
23 
24   You may use `SNESVISetVariableBounds()` to provide the bounds once if they will never change
25 
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.
28 
29 .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `DMSetVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`,
30           `'SNESSetType()`
31 @*/
32 PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES snes, Vec lower, Vec higher))
33 {
34   PetscErrorCode (*f)(SNES, PetscErrorCode (*)(SNES, Vec, Vec));
35 
36   PetscFunctionBegin;
37   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
38   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", &f));
39   if (f) PetscUseMethod(snes, "SNESVISetComputeVariableBounds_C", (SNES, PetscErrorCode(*)(SNES, Vec, Vec)), (snes, compute));
40   else PetscCall(SNESVISetComputeVariableBounds_VI(snes, compute));
41   PetscFunctionReturn(PETSC_SUCCESS);
42 }
43 
44 PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes, SNESVIComputeVariableBoundsFunction compute)
45 {
46   PetscFunctionBegin;
47   snes->ops->computevariablebounds = compute;
48   PetscFunctionReturn(PETSC_SUCCESS);
49 }
50 
51 PetscErrorCode SNESVIMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy)
52 {
53   Vec         X, F, Finactive;
54   IS          isactive;
55   PetscViewer viewer = (PetscViewer)dummy;
56 
57   PetscFunctionBegin;
58   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
59   PetscCall(SNESGetSolution(snes, &X));
60   PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive));
61   PetscCall(VecDuplicate(F, &Finactive));
62   PetscCall(VecCopy(F, Finactive));
63   PetscCall(VecISSet(Finactive, isactive, 0.0));
64   PetscCall(ISDestroy(&isactive));
65   PetscCall(VecView(Finactive, viewer));
66   PetscCall(VecDestroy(&Finactive));
67   PetscFunctionReturn(PETSC_SUCCESS);
68 }
69 
70 PetscErrorCode SNESMonitorVI(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy)
71 {
72   PetscViewer        viewer = (PetscViewer)dummy;
73   const PetscScalar *x, *xl, *xu, *f;
74   PetscInt           i, n, act[2] = {0, 0}, fact[2], N;
75   /* Number of components that actually hit the bounds (c.f. active variables) */
76   PetscInt  act_bound[2] = {0, 0}, fact_bound[2];
77   PetscReal rnorm, fnorm, zerotolerance = snes->vizerotolerance;
78   double    tmp;
79 
80   PetscFunctionBegin;
81   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
82   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
83   PetscCall(VecGetSize(snes->vec_sol, &N));
84   PetscCall(VecGetArrayRead(snes->xl, &xl));
85   PetscCall(VecGetArrayRead(snes->xu, &xu));
86   PetscCall(VecGetArrayRead(snes->vec_sol, &x));
87   PetscCall(VecGetArrayRead(snes->vec_func, &f));
88 
89   rnorm = 0.0;
90   for (i = 0; i < n; i++) {
91     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]);
92     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++;
93     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++;
94     else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Can never get here");
95   }
96 
97   for (i = 0; i < n; i++) {
98     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++;
99     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++;
100   }
101   PetscCall(VecRestoreArrayRead(snes->vec_func, &f));
102   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
103   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
104   PetscCall(VecRestoreArrayRead(snes->vec_sol, &x));
105   PetscCall(MPIU_Allreduce(&rnorm, &fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
106   PetscCall(MPIU_Allreduce(act, fact, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
107   PetscCall(MPIU_Allreduce(act_bound, fact_bound, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
108   fnorm = PetscSqrtReal(fnorm);
109 
110   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
111   if (snes->ntruebounds) tmp = ((double)(fact[0] + fact[1])) / ((double)snes->ntruebounds);
112   else tmp = 0.0;
113   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));
114 
115   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
116   PetscFunctionReturn(PETSC_SUCCESS);
117 }
118 
119 /*
120      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
121     || 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
122     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
123     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
124 */
125 PetscErrorCode SNESVICheckLocalMin_Private(SNES snes, Mat A, Vec F, Vec W, PetscReal fnorm, PetscBool *ismin)
126 {
127   PetscReal a1;
128   PetscBool hastranspose;
129 
130   PetscFunctionBegin;
131   *ismin = PETSC_FALSE;
132   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
133   if (hastranspose) {
134     /* Compute || J^T F|| */
135     PetscCall(MatMultTranspose(A, F, W));
136     PetscCall(VecNorm(W, NORM_2, &a1));
137     PetscCall(PetscInfo(snes, "|| J^T F|| %g near zero implies found a local minimum\n", (double)(a1 / fnorm)));
138     if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE;
139   } else {
140     Vec         work;
141     PetscScalar result;
142     PetscReal   wnorm;
143 
144     PetscCall(VecSetRandom(W, NULL));
145     PetscCall(VecNorm(W, NORM_2, &wnorm));
146     PetscCall(VecDuplicate(W, &work));
147     PetscCall(MatMult(A, W, work));
148     PetscCall(VecDot(F, work, &result));
149     PetscCall(VecDestroy(&work));
150     a1 = PetscAbsScalar(result) / (fnorm * wnorm);
151     PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n", (double)a1));
152     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
153   }
154   PetscFunctionReturn(PETSC_SUCCESS);
155 }
156 
157 /*
158      Checks if J^T(F - J*X) = 0
159 */
160 PetscErrorCode SNESVICheckResidual_Private(SNES snes, Mat A, Vec F, Vec X, Vec W1, Vec W2)
161 {
162   PetscReal a1, a2;
163   PetscBool hastranspose;
164 
165   PetscFunctionBegin;
166   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
167   if (hastranspose) {
168     PetscCall(MatMult(A, X, W1));
169     PetscCall(VecAXPY(W1, -1.0, F));
170 
171     /* Compute || J^T W|| */
172     PetscCall(MatMultTranspose(A, W1, W2));
173     PetscCall(VecNorm(W1, NORM_2, &a1));
174     PetscCall(VecNorm(W2, NORM_2, &a2));
175     if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n", (double)(a2 / a1)));
176   }
177   PetscFunctionReturn(PETSC_SUCCESS);
178 }
179 
180 /*
181   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
182 
183   Notes:
184   The convergence criterion currently implemented is
185   merit < abstol
186   merit < rtol*merit_initial
187 */
188 PetscErrorCode SNESConvergedDefault_VI(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gradnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
189 {
190   PetscFunctionBegin;
191   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
192   PetscAssertPointer(reason, 6);
193 
194   *reason = SNES_CONVERGED_ITERATING;
195 
196   if (!it) {
197     /* set parameter for default relative tolerance convergence test */
198     snes->ttol = fnorm * snes->rtol;
199   }
200   if (fnorm != fnorm) {
201     PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n"));
202     *reason = SNES_DIVERGED_FNORM_NAN;
203   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
204     PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g\n", (double)fnorm, (double)snes->abstol));
205     *reason = SNES_CONVERGED_FNORM_ABS;
206   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
207     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs));
208     *reason = SNES_DIVERGED_FUNCTION_COUNT;
209   }
210 
211   if (it && !*reason) {
212     if (fnorm < snes->ttol) {
213       PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g (relative tolerance)\n", (double)fnorm, (double)snes->ttol));
214       *reason = SNES_CONVERGED_FNORM_RELATIVE;
215     }
216   }
217   PetscFunctionReturn(PETSC_SUCCESS);
218 }
219 
220 /*
221    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
222 
223    Input Parameters:
224 .  SNES - nonlinear solver context
225 
226    Output Parameters:
227 .  X - Bound projected X
228 
229 */
230 
231 PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X)
232 {
233   const PetscScalar *xl, *xu;
234   PetscScalar       *x;
235   PetscInt           i, n;
236 
237   PetscFunctionBegin;
238   PetscCall(VecGetLocalSize(X, &n));
239   PetscCall(VecGetArray(X, &x));
240   PetscCall(VecGetArrayRead(snes->xl, &xl));
241   PetscCall(VecGetArrayRead(snes->xu, &xu));
242 
243   for (i = 0; i < n; i++) {
244     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
245     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
246   }
247   PetscCall(VecRestoreArray(X, &x));
248   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
249   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
250   PetscFunctionReturn(PETSC_SUCCESS);
251 }
252 
253 /*@
254   SNESVIGetActiveSetIS - Gets the global indices for the active set variables
255 
256   Input Parameters:
257 + snes - the `SNES` context
258 . X    - the `snes` solution vector
259 - F    - the nonlinear function vector
260 
261   Output Parameter:
262 . ISact - active set index set
263 
264   Level: developer
265 
266 .seealso: `SNES`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`
267 @*/
268 PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact)
269 {
270   Vec                Xl = snes->xl, Xu = snes->xu;
271   const PetscScalar *x, *f, *xl, *xu;
272   PetscInt          *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0;
273   PetscReal          zerotolerance = snes->vizerotolerance;
274 
275   PetscFunctionBegin;
276   PetscCall(VecGetLocalSize(X, &nlocal));
277   PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh));
278   PetscCall(VecGetArrayRead(X, &x));
279   PetscCall(VecGetArrayRead(Xl, &xl));
280   PetscCall(VecGetArrayRead(Xu, &xu));
281   PetscCall(VecGetArrayRead(F, &f));
282   /* Compute active set size */
283   for (i = 0; i < nlocal; i++) {
284     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++;
285   }
286 
287   PetscCall(PetscMalloc1(nloc_isact, &idx_act));
288 
289   /* Set active set indices */
290   for (i = 0; i < nlocal; i++) {
291     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;
292   }
293 
294   /* Create active set IS */
295   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact));
296 
297   PetscCall(VecRestoreArrayRead(X, &x));
298   PetscCall(VecRestoreArrayRead(Xl, &xl));
299   PetscCall(VecRestoreArrayRead(Xu, &xu));
300   PetscCall(VecRestoreArrayRead(F, &f));
301   PetscFunctionReturn(PETSC_SUCCESS);
302 }
303 
304 PetscErrorCode SNESVICreateIndexSets_RS(SNES snes, Vec X, Vec F, IS *ISact, IS *ISinact)
305 {
306   PetscInt rstart, rend;
307 
308   PetscFunctionBegin;
309   PetscCall(SNESVIGetActiveSetIS(snes, X, F, ISact));
310   PetscCall(VecGetOwnershipRange(X, &rstart, &rend));
311   PetscCall(ISComplement(*ISact, rstart, rend, ISinact));
312   PetscFunctionReturn(PETSC_SUCCESS);
313 }
314 
315 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm)
316 {
317   const PetscScalar *x, *xl, *xu, *f;
318   PetscInt           i, n;
319   PetscReal          rnorm, zerotolerance = snes->vizerotolerance;
320 
321   PetscFunctionBegin;
322   PetscCall(VecGetLocalSize(X, &n));
323   PetscCall(VecGetArrayRead(snes->xl, &xl));
324   PetscCall(VecGetArrayRead(snes->xu, &xu));
325   PetscCall(VecGetArrayRead(X, &x));
326   PetscCall(VecGetArrayRead(F, &f));
327   rnorm = 0.0;
328   for (i = 0; i < n; i++) {
329     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]);
330   }
331   PetscCall(VecRestoreArrayRead(F, &f));
332   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
333   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
334   PetscCall(VecRestoreArrayRead(X, &x));
335   PetscCall(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
336   *fnorm = PetscSqrtReal(*fnorm);
337   PetscFunctionReturn(PETSC_SUCCESS);
338 }
339 
340 PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu)
341 {
342   PetscFunctionBegin;
343   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
344   PetscFunctionReturn(PETSC_SUCCESS);
345 }
346 
347 /*
348    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
349    of the SNESVI nonlinear solver.
350 
351    Input Parameter:
352 .  snes - the SNES context
353 
354    Application Interface Routine: SNESSetUp()
355 
356    Notes:
357    For basic use of the SNES solvers, the user need not explicitly call
358    SNESSetUp(), since these actions will automatically occur during
359    the call to SNESSolve().
360  */
361 PetscErrorCode SNESSetUp_VI(SNES snes)
362 {
363   PetscInt i_start[3], i_end[3];
364 
365   PetscFunctionBegin;
366   PetscCall(SNESSetWorkVecs(snes, 1));
367   PetscCall(SNESSetUpMatrices(snes));
368 
369   if (!snes->ops->computevariablebounds && snes->dm) {
370     PetscBool flag;
371     PetscCall(DMHasVariableBounds(snes->dm, &flag));
372     if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
373   }
374   if (!snes->usersetbounds) {
375     if (snes->ops->computevariablebounds) {
376       if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
377       if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
378       PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
379     } else if (!snes->xl && !snes->xu) {
380       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
381       PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
382       PetscCall(VecSet(snes->xl, PETSC_NINFINITY));
383       PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
384       PetscCall(VecSet(snes->xu, PETSC_INFINITY));
385     } else {
386       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
387       PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end));
388       PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1));
389       PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2));
390       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]))
391         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.");
392     }
393   }
394   PetscFunctionReturn(PETSC_SUCCESS);
395 }
396 PetscErrorCode SNESReset_VI(SNES snes)
397 {
398   PetscFunctionBegin;
399   PetscCall(VecDestroy(&snes->xl));
400   PetscCall(VecDestroy(&snes->xu));
401   snes->usersetbounds = PETSC_FALSE;
402   PetscFunctionReturn(PETSC_SUCCESS);
403 }
404 
405 /*
406    SNESDestroy_VI - Destroys the private SNES_VI context that was created
407    with SNESCreate_VI().
408 
409    Input Parameter:
410 .  snes - the SNES context
411 
412    Application Interface Routine: SNESDestroy()
413  */
414 PetscErrorCode SNESDestroy_VI(SNES snes)
415 {
416   PetscFunctionBegin;
417   PetscCall(PetscFree(snes->data));
418 
419   /* clear composed functions */
420   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL));
421   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL));
422   PetscFunctionReturn(PETSC_SUCCESS);
423 }
424 
425 /*@
426   SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving
427   (differential) variable inequalities.
428 
429   Input Parameters:
430 + snes - the `SNES` context.
431 . xl   - lower bound.
432 - xu   - upper bound.
433 
434   Level: advanced
435 
436   Notes:
437   If this routine is not called then the lower and upper bounds are set to
438   `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`.
439 
440   Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.
441 
442   For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY`
443 
444   `SNESVISetComputeVariableBounds()` can be used to provide a function that computes the bounds. This should be used if you are using, for example, grid
445   sequencing and need bounds set for a variety of vectors
446 
447 .seealso: [](sec_vi), `SNES`, `SNESVIGetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `'SNESSetType()`
448 @*/
449 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
450 {
451   PetscErrorCode (*f)(SNES, Vec, Vec);
452 
453   PetscFunctionBegin;
454   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
455   PetscValidHeaderSpecific(xl, VEC_CLASSID, 2);
456   PetscValidHeaderSpecific(xu, VEC_CLASSID, 3);
457   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f));
458   if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu));
459   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
460   snes->usersetbounds = PETSC_TRUE;
461   PetscFunctionReturn(PETSC_SUCCESS);
462 }
463 
464 PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu)
465 {
466   const PetscScalar *xxl, *xxu;
467   PetscInt           i, n, cnt = 0;
468 
469   PetscFunctionBegin;
470   PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL));
471   PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first");
472   {
473     PetscInt xlN, xuN, N;
474     PetscCall(VecGetSize(xl, &xlN));
475     PetscCall(VecGetSize(xu, &xuN));
476     PetscCall(VecGetSize(snes->vec_func, &N));
477     PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N);
478     PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N);
479   }
480   PetscCall(PetscObjectReference((PetscObject)xl));
481   PetscCall(PetscObjectReference((PetscObject)xu));
482   PetscCall(VecDestroy(&snes->xl));
483   PetscCall(VecDestroy(&snes->xu));
484   snes->xl = xl;
485   snes->xu = xu;
486   PetscCall(VecGetLocalSize(xl, &n));
487   PetscCall(VecGetArrayRead(xl, &xxl));
488   PetscCall(VecGetArrayRead(xu, &xxu));
489   for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
490 
491   PetscCall(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
492   PetscCall(VecRestoreArrayRead(xl, &xxl));
493   PetscCall(VecRestoreArrayRead(xu, &xxu));
494   PetscFunctionReturn(PETSC_SUCCESS);
495 }
496 
497 /*@
498   SNESVIGetVariableBounds - Gets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving
499   (differential) variable inequalities.
500 
501   Input Parameters:
502 + snes - the `SNES` context.
503 . xl   - lower bound (may be `NULL`)
504 - xu   - upper bound (may be `NULL`)
505 
506   Level: advanced
507 
508   Notes:
509   These vectors are owned by the `SNESVI` and should not be destroyed by the caller
510 
511 .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `'SNESSetType()`
512 @*/
513 PetscErrorCode SNESVIGetVariableBounds(SNES snes, Vec *xl, Vec *xu)
514 {
515   PetscFunctionBegin;
516   PetscCheck(snes->usersetbounds, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must set SNESVI bounds before calling SNESVIGetVariableBounds()");
517   if (xl) *xl = snes->xl;
518   if (xu) *xu = snes->xu;
519   PetscFunctionReturn(PETSC_SUCCESS);
520 }
521 
522 PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject)
523 {
524   PetscBool flg = PETSC_FALSE;
525 
526   PetscFunctionBegin;
527   PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options");
528   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL));
529   PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL));
530   if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL));
531   flg = PETSC_FALSE;
532   PetscCall(PetscOptionsBool("-snes_vi_monitor_residual", "Monitor residual all non-active variables; using zero for active constraints", "SNESMonitorVIResidual", flg, &flg, NULL));
533   if (flg) PetscCall(SNESMonitorSet(snes, SNESVIMonitorResidual, PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)), NULL));
534   PetscOptionsHeadEnd();
535   PetscFunctionReturn(PETSC_SUCCESS);
536 }
537