xref: /petsc/src/snes/impls/vi/vi.c (revision bcd3bd92eda2d5998e2f14c4bbfb33bd936bdc3e)
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 Parameter:
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 PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact)
265 {
266   Vec                Xl = snes->xl, Xu = snes->xu;
267   const PetscScalar *x, *f, *xl, *xu;
268   PetscInt          *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0;
269   PetscReal          zerotolerance = snes->vizerotolerance;
270 
271   PetscFunctionBegin;
272   PetscCall(VecGetLocalSize(X, &nlocal));
273   PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh));
274   PetscCall(VecGetArrayRead(X, &x));
275   PetscCall(VecGetArrayRead(Xl, &xl));
276   PetscCall(VecGetArrayRead(Xu, &xu));
277   PetscCall(VecGetArrayRead(F, &f));
278   /* Compute active set size */
279   for (i = 0; i < nlocal; i++) {
280     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++;
281   }
282 
283   PetscCall(PetscMalloc1(nloc_isact, &idx_act));
284 
285   /* Set active set indices */
286   for (i = 0; i < nlocal; i++) {
287     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;
288   }
289 
290   /* Create active set IS */
291   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact));
292 
293   PetscCall(VecRestoreArrayRead(X, &x));
294   PetscCall(VecRestoreArrayRead(Xl, &xl));
295   PetscCall(VecRestoreArrayRead(Xu, &xu));
296   PetscCall(VecRestoreArrayRead(F, &f));
297   PetscFunctionReturn(PETSC_SUCCESS);
298 }
299 
300 PetscErrorCode SNESVICreateIndexSets_RS(SNES snes, Vec X, Vec F, IS *ISact, IS *ISinact)
301 {
302   PetscInt rstart, rend;
303 
304   PetscFunctionBegin;
305   PetscCall(SNESVIGetActiveSetIS(snes, X, F, ISact));
306   PetscCall(VecGetOwnershipRange(X, &rstart, &rend));
307   PetscCall(ISComplement(*ISact, rstart, rend, ISinact));
308   PetscFunctionReturn(PETSC_SUCCESS);
309 }
310 
311 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm)
312 {
313   const PetscScalar *x, *xl, *xu, *f;
314   PetscInt           i, n;
315   PetscReal          rnorm, zerotolerance = snes->vizerotolerance;
316 
317   PetscFunctionBegin;
318   PetscCall(VecGetLocalSize(X, &n));
319   PetscCall(VecGetArrayRead(snes->xl, &xl));
320   PetscCall(VecGetArrayRead(snes->xu, &xu));
321   PetscCall(VecGetArrayRead(X, &x));
322   PetscCall(VecGetArrayRead(F, &f));
323   rnorm = 0.0;
324   for (i = 0; i < n; i++) {
325     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]);
326   }
327   PetscCall(VecRestoreArrayRead(F, &f));
328   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
329   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
330   PetscCall(VecRestoreArrayRead(X, &x));
331   PetscCall(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
332   *fnorm = PetscSqrtReal(*fnorm);
333   PetscFunctionReturn(PETSC_SUCCESS);
334 }
335 
336 PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu)
337 {
338   PetscFunctionBegin;
339   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
340   PetscFunctionReturn(PETSC_SUCCESS);
341 }
342 
343 /*
344    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
345    of the SNESVI nonlinear solver.
346 
347    Input Parameter:
348 .  snes - the SNES context
349 
350    Application Interface Routine: SNESSetUp()
351 
352    Notes:
353    For basic use of the SNES solvers, the user need not explicitly call
354    SNESSetUp(), since these actions will automatically occur during
355    the call to SNESSolve().
356  */
357 PetscErrorCode SNESSetUp_VI(SNES snes)
358 {
359   PetscInt i_start[3], i_end[3];
360 
361   PetscFunctionBegin;
362   PetscCall(SNESSetWorkVecs(snes, 1));
363   PetscCall(SNESSetUpMatrices(snes));
364 
365   if (!snes->ops->computevariablebounds && snes->dm) {
366     PetscBool flag;
367     PetscCall(DMHasVariableBounds(snes->dm, &flag));
368     if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
369   }
370   if (!snes->usersetbounds) {
371     if (snes->ops->computevariablebounds) {
372       if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
373       if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
374       PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
375     } else if (!snes->xl && !snes->xu) {
376       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
377       PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
378       PetscCall(VecSet(snes->xl, PETSC_NINFINITY));
379       PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
380       PetscCall(VecSet(snes->xu, PETSC_INFINITY));
381     } else {
382       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
383       PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end));
384       PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1));
385       PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2));
386       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]))
387         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.");
388     }
389   }
390   PetscFunctionReturn(PETSC_SUCCESS);
391 }
392 PetscErrorCode SNESReset_VI(SNES snes)
393 {
394   PetscFunctionBegin;
395   PetscCall(VecDestroy(&snes->xl));
396   PetscCall(VecDestroy(&snes->xu));
397   snes->usersetbounds = PETSC_FALSE;
398   PetscFunctionReturn(PETSC_SUCCESS);
399 }
400 
401 /*
402    SNESDestroy_VI - Destroys the private SNES_VI context that was created
403    with SNESCreate_VI().
404 
405    Input Parameter:
406 .  snes - the SNES context
407 
408    Application Interface Routine: SNESDestroy()
409  */
410 PetscErrorCode SNESDestroy_VI(SNES snes)
411 {
412   PetscFunctionBegin;
413   PetscCall(PetscFree(snes->data));
414 
415   /* clear composed functions */
416   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL));
417   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL));
418   PetscFunctionReturn(PETSC_SUCCESS);
419 }
420 
421 /*@
422   SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving
423   (differential) variable inequalities.
424 
425   Input Parameters:
426 + snes - the `SNES` context.
427 . xl   - lower bound.
428 - xu   - upper bound.
429 
430   Level: advanced
431 
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()`.
435 
436   Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.
437 
438   For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY`
439 
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
442 
443 .seealso: [](sec_vi), `SNES`, `SNESVIGetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `'SNESSetType()`
444 @*/
445 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
446 {
447   PetscErrorCode (*f)(SNES, Vec, Vec);
448 
449   PetscFunctionBegin;
450   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
451   PetscValidHeaderSpecific(xl, VEC_CLASSID, 2);
452   PetscValidHeaderSpecific(xu, VEC_CLASSID, 3);
453   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f));
454   if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu));
455   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
456   snes->usersetbounds = PETSC_TRUE;
457   PetscFunctionReturn(PETSC_SUCCESS);
458 }
459 
460 PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu)
461 {
462   const PetscScalar *xxl, *xxu;
463   PetscInt           i, n, cnt = 0;
464 
465   PetscFunctionBegin;
466   PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL));
467   PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first");
468   {
469     PetscInt xlN, xuN, N;
470     PetscCall(VecGetSize(xl, &xlN));
471     PetscCall(VecGetSize(xu, &xuN));
472     PetscCall(VecGetSize(snes->vec_func, &N));
473     PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N);
474     PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N);
475   }
476   PetscCall(PetscObjectReference((PetscObject)xl));
477   PetscCall(PetscObjectReference((PetscObject)xu));
478   PetscCall(VecDestroy(&snes->xl));
479   PetscCall(VecDestroy(&snes->xu));
480   snes->xl = xl;
481   snes->xu = xu;
482   PetscCall(VecGetLocalSize(xl, &n));
483   PetscCall(VecGetArrayRead(xl, &xxl));
484   PetscCall(VecGetArrayRead(xu, &xxu));
485   for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
486 
487   PetscCall(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
488   PetscCall(VecRestoreArrayRead(xl, &xxl));
489   PetscCall(VecRestoreArrayRead(xu, &xxu));
490   PetscFunctionReturn(PETSC_SUCCESS);
491 }
492 
493 /*@
494   SNESVIGetVariableBounds - Gets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving
495   (differential) variable inequalities.
496 
497   Input Parameters:
498 + snes - the `SNES` context.
499 . xl   - lower bound (may be `NULL`)
500 - xu   - upper bound (may be `NULL`)
501 
502   Level: advanced
503 
504   Notes:
505   These vectors are owned by the `SNESVI` and should not be destroyed by the caller
506 
507 .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `'SNESSetType()`
508 @*/
509 PetscErrorCode SNESVIGetVariableBounds(SNES snes, Vec *xl, Vec *xu)
510 {
511   PetscFunctionBegin;
512   PetscCheck(snes->usersetbounds, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must set SNESVI bounds before calling SNESVIGetVariableBounds()");
513   if (xl) *xl = snes->xl;
514   if (xu) *xu = snes->xu;
515   PetscFunctionReturn(PETSC_SUCCESS);
516 }
517 
518 PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject)
519 {
520   PetscBool flg = PETSC_FALSE;
521 
522   PetscFunctionBegin;
523   PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options");
524   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL));
525   PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL));
526   if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL));
527   flg = PETSC_FALSE;
528   PetscCall(PetscOptionsBool("-snes_vi_monitor_residual", "Monitor residual all non-active variables; using zero for active constraints", "SNESMonitorVIResidual", flg, &flg, NULL));
529   if (flg) PetscCall(SNESMonitorSet(snes, SNESVIMonitorResidual, PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)), NULL));
530   PetscOptionsHeadEnd();
531   PetscFunctionReturn(PETSC_SUCCESS);
532 }
533