xref: /petsc/src/snes/impls/vi/vi.c (revision bd89dbf26d8a5efecb980364933175da61864cd7)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
21e25c274SJed Brown #include <petscdm.h>
3d0af7cd3SBarry Smith 
42176524cSBarry Smith /*@C
5f6dfbefdSBarry Smith   SNESVISetComputeVariableBounds - Sets a function that is called to compute the bounds on variable for
6f6dfbefdSBarry Smith   (differential) variable inequalities.
72176524cSBarry Smith 
8e4094ef1SJacob Faibussowitsch   Input Parameters:
9f6dfbefdSBarry Smith + snes    - the `SNES` context
10f6dfbefdSBarry Smith - compute - function that computes the bounds
11f6dfbefdSBarry Smith 
12e4094ef1SJacob Faibussowitsch   Calling sequence of `compute`:
13f6dfbefdSBarry Smith + snes   - the `SNES` context
14f6dfbefdSBarry Smith . lower  - vector to hold lower bounds
15f6dfbefdSBarry Smith - higher - vector to hold upper bounds
162176524cSBarry Smith 
172bd2b0e6SSatish Balay   Level: advanced
182176524cSBarry Smith 
19f6dfbefdSBarry Smith   Notes:
20f0b84518SBarry Smith   Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.
21f0b84518SBarry Smith 
22f6dfbefdSBarry Smith   For entries with no bounds you can set `PETSC_NINFINITY` or `PETSC_INFINITY`
23c1c3a0ecSBarry Smith 
24f6dfbefdSBarry Smith   You may use `SNESVISetVariableBounds()` to provide the bounds once if they will never change
25f6dfbefdSBarry Smith 
26f6dfbefdSBarry Smith   If you have associated a `DM` with the `SNES` and provided a function to the `DM` via `DMSetVariableBounds()` that will be used automatically
27f6dfbefdSBarry Smith   to provide the bounds and you need not use this function.
28f6dfbefdSBarry Smith 
2933e0caf2SBarry Smith .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `DMSetVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`,
302913281cSPierre Jolivet           `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
31aab9d709SJed Brown @*/
SNESVISetComputeVariableBounds(SNES snes,PetscErrorCode (* compute)(SNES snes,Vec lower,Vec higher))32e4094ef1SJacob Faibussowitsch PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES snes, Vec lower, Vec higher))
33d71ae5a4SJacob Faibussowitsch {
345f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES, PetscErrorCode (*)(SNES, Vec, Vec));
352176524cSBarry Smith 
362176524cSBarry Smith   PetscFunctionBegin;
3761589011SJed Brown   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
389566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", &f));
39cac4c232SBarry Smith   if (f) PetscUseMethod(snes, "SNESVISetComputeVariableBounds_C", (SNES, PetscErrorCode (*)(SNES, Vec, Vec)), (snes, compute));
409566063dSJacob Faibussowitsch   else PetscCall(SNESVISetComputeVariableBounds_VI(snes, compute));
413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
422176524cSBarry Smith }
432176524cSBarry Smith 
SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFn * compute)448434afd1SBarry Smith PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes, SNESVIComputeVariableBoundsFn *compute)
45d71ae5a4SJacob Faibussowitsch {
4661589011SJed Brown   PetscFunctionBegin;
4761589011SJed Brown   snes->ops->computevariablebounds = compute;
483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4961589011SJed Brown }
502176524cSBarry Smith 
SNESVIMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat * vf)51b6f515dbSMatthew G. Knepley static PetscErrorCode SNESVIMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
52d71ae5a4SJacob Faibussowitsch {
53ffdf2a20SBarry Smith   Vec X, F, Finactive;
54ffdf2a20SBarry Smith   IS  isactive;
55ffdf2a20SBarry Smith 
56ffdf2a20SBarry Smith   PetscFunctionBegin;
57b6f515dbSMatthew G. Knepley   PetscValidHeaderSpecific(vf->viewer, PETSC_VIEWER_CLASSID, 4);
589566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
599566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &X));
609566063dSJacob Faibussowitsch   PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive));
619566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(F, &Finactive));
62b6f515dbSMatthew G. Knepley   PetscCall(PetscObjectCompose((PetscObject)Finactive, "__Vec_bc_zero__", (PetscObject)snes));
639566063dSJacob Faibussowitsch   PetscCall(VecCopy(F, Finactive));
649566063dSJacob Faibussowitsch   PetscCall(VecISSet(Finactive, isactive, 0.0));
659566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isactive));
66b6f515dbSMatthew G. Knepley   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
67b6f515dbSMatthew G. Knepley   PetscCall(VecView(Finactive, vf->viewer));
68b6f515dbSMatthew G. Knepley   PetscCall(PetscViewerPopFormat(vf->viewer));
69b6f515dbSMatthew G. Knepley   PetscCall(PetscObjectCompose((PetscObject)Finactive, "__Vec_bc_zero__", NULL));
709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Finactive));
713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72ffdf2a20SBarry Smith }
73ffdf2a20SBarry Smith 
SNESVIMonitorActive(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat * vf)74b6f515dbSMatthew G. Knepley static PetscErrorCode SNESVIMonitorActive(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
75b6f515dbSMatthew G. Knepley {
76b6f515dbSMatthew G. Knepley   Vec X, F, A;
77b6f515dbSMatthew G. Knepley   IS  isactive;
78b6f515dbSMatthew G. Knepley 
79b6f515dbSMatthew G. Knepley   PetscFunctionBegin;
80b6f515dbSMatthew G. Knepley   PetscValidHeaderSpecific(vf->viewer, PETSC_VIEWER_CLASSID, 4);
81b6f515dbSMatthew G. Knepley   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
82b6f515dbSMatthew G. Knepley   PetscCall(SNESGetSolution(snes, &X));
83b6f515dbSMatthew G. Knepley   PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive));
84b6f515dbSMatthew G. Knepley   PetscCall(VecDuplicate(F, &A));
85b6f515dbSMatthew G. Knepley   PetscCall(PetscObjectCompose((PetscObject)A, "__Vec_bc_zero__", (PetscObject)snes));
86b6f515dbSMatthew G. Knepley   PetscCall(VecSet(A, 0.));
87b6f515dbSMatthew G. Knepley   PetscCall(VecISSet(A, isactive, 1.));
88b6f515dbSMatthew G. Knepley   PetscCall(ISDestroy(&isactive));
89b6f515dbSMatthew G. Knepley   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
90b6f515dbSMatthew G. Knepley   PetscCall(VecView(A, vf->viewer));
91b6f515dbSMatthew G. Knepley   PetscCall(PetscViewerPopFormat(vf->viewer));
92b6f515dbSMatthew G. Knepley   PetscCall(PetscObjectCompose((PetscObject)A, "__Vec_bc_zero__", NULL));
93b6f515dbSMatthew G. Knepley   PetscCall(VecDestroy(&A));
94b6f515dbSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
95b6f515dbSMatthew G. Knepley }
96b6f515dbSMatthew G. Knepley 
SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void * dummy)97ba38deedSJacob Faibussowitsch static PetscErrorCode SNESMonitorVI(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy)
98d71ae5a4SJacob Faibussowitsch {
994d4332d5SBarry Smith   PetscViewer        viewer = (PetscViewer)dummy;
1009308c137SBarry Smith   const PetscScalar *x, *xl, *xu, *f;
1016fd67740SBarry Smith   PetscInt           i, n, act[2] = {0, 0}, fact[2], N;
1026a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
1036a9e2d46SJungho Lee   PetscInt  act_bound[2] = {0, 0}, fact_bound[2];
104349187a7SBarry Smith   PetscReal rnorm, fnorm, zerotolerance = snes->vizerotolerance;
1059d1809e2SSatish Balay   double    tmp;
1069308c137SBarry Smith 
1079308c137SBarry Smith   PetscFunctionBegin;
1084d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
1099566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
1109566063dSJacob Faibussowitsch   PetscCall(VecGetSize(snes->vec_sol, &N));
1119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
1129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
1139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->vec_sol, &x));
1149566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->vec_func, &f));
1159308c137SBarry Smith 
1169308c137SBarry Smith   rnorm = 0.0;
1179308c137SBarry Smith   for (i = 0; i < n; i++) {
118f4f49eeaSPierre Jolivet     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]);
119349187a7SBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++;
120349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++;
121ce94432eSBarry Smith     else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Can never get here");
1229308c137SBarry Smith   }
1236a9e2d46SJungho Lee 
1246a9e2d46SJungho Lee   for (i = 0; i < n; i++) {
125349187a7SBarry Smith     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++;
126349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++;
1276a9e2d46SJungho Lee   }
1289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->vec_func, &f));
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
1309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
1319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->vec_sol, &x));
132462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&rnorm, &fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
133462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(act, fact, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
134462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(act_bound, fact_bound, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
135f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1366fd67740SBarry Smith 
1379566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
1389d1809e2SSatish Balay   if (snes->ntruebounds) tmp = ((double)(fact[0] + fact[1])) / ((double)snes->ntruebounds);
1399d1809e2SSatish Balay   else tmp = 0.0;
14063a3b9bcSJacob Faibussowitsch   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));
1416a9e2d46SJungho Lee 
1429566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
1433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1449308c137SBarry Smith }
1459308c137SBarry Smith 
1462b4ed282SShri Abhyankar /*
1472b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1482b4ed282SShri Abhyankar     || 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
1492b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1502b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1512b4ed282SShri Abhyankar */
SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool * ismin)152d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVICheckLocalMin_Private(SNES snes, Mat A, Vec F, Vec W, PetscReal fnorm, PetscBool *ismin)
153d71ae5a4SJacob Faibussowitsch {
1542b4ed282SShri Abhyankar   PetscReal a1;
155ace3abfcSBarry Smith   PetscBool hastranspose;
1562b4ed282SShri Abhyankar 
1572b4ed282SShri Abhyankar   PetscFunctionBegin;
1582b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1599566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
1602b4ed282SShri Abhyankar   if (hastranspose) {
1612b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1629566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, F, W));
1639566063dSJacob Faibussowitsch     PetscCall(VecNorm(W, NORM_2, &a1));
1649566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "|| J^T F|| %g near zero implies found a local minimum\n", (double)(a1 / fnorm)));
1652b4ed282SShri Abhyankar     if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE;
1662b4ed282SShri Abhyankar   } else {
1672b4ed282SShri Abhyankar     Vec         work;
1682b4ed282SShri Abhyankar     PetscScalar result;
1692b4ed282SShri Abhyankar     PetscReal   wnorm;
1702b4ed282SShri Abhyankar 
1719566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(W, NULL));
1729566063dSJacob Faibussowitsch     PetscCall(VecNorm(W, NORM_2, &wnorm));
1739566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(W, &work));
1749566063dSJacob Faibussowitsch     PetscCall(MatMult(A, W, work));
1759566063dSJacob Faibussowitsch     PetscCall(VecDot(F, work, &result));
1769566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&work));
1772b4ed282SShri Abhyankar     a1 = PetscAbsScalar(result) / (fnorm * wnorm);
1789566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n", (double)a1));
1792b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1802b4ed282SShri Abhyankar   }
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1822b4ed282SShri Abhyankar }
1832b4ed282SShri Abhyankar 
1842b4ed282SShri Abhyankar /*
1858d359177SBarry Smith   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
1862b4ed282SShri Abhyankar 
1872b4ed282SShri Abhyankar   Notes:
1882b4ed282SShri Abhyankar   The convergence criterion currently implemented is
1892b4ed282SShri Abhyankar   merit < abstol
1902b4ed282SShri Abhyankar   merit < rtol*merit_initial
1912b4ed282SShri Abhyankar */
SNESConvergedDefault_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason * reason,void * dummy)192d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESConvergedDefault_VI(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gradnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
193d71ae5a4SJacob Faibussowitsch {
1942b4ed282SShri Abhyankar   PetscFunctionBegin;
1952b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1964f572ea9SToby Isaac   PetscAssertPointer(reason, 6);
1972b4ed282SShri Abhyankar 
1982b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
1992b4ed282SShri Abhyankar 
2002b4ed282SShri Abhyankar   if (!it) {
2012b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
2027fe79bc4SShri Abhyankar     snes->ttol = fnorm * snes->rtol;
2032b4ed282SShri Abhyankar   }
2047fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
2059566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n"));
206*76c63389SBarry Smith     *reason = SNES_DIVERGED_FUNCTION_NANORINF;
2070d6f27a8SBarry Smith   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
2089566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g\n", (double)fnorm, (double)snes->abstol));
2092b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
210e71169deSBarry Smith   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
21163a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs));
2122b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
2132b4ed282SShri Abhyankar   }
2142b4ed282SShri Abhyankar 
2152b4ed282SShri Abhyankar   if (it && !*reason) {
2167fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2179566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g (relative tolerance)\n", (double)fnorm, (double)snes->ttol));
2182b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2192b4ed282SShri Abhyankar     }
2202b4ed282SShri Abhyankar   }
2213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2222b4ed282SShri Abhyankar }
2232b4ed282SShri Abhyankar 
224c1a5e521SShri Abhyankar /*
225c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
226c1a5e521SShri Abhyankar 
227c1a5e521SShri Abhyankar    Input Parameters:
228c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
229c1a5e521SShri Abhyankar 
230c1a5e521SShri Abhyankar    Output Parameters:
231c1a5e521SShri Abhyankar .  X - Bound projected X
232c1a5e521SShri Abhyankar 
233c1a5e521SShri Abhyankar */
234c1a5e521SShri Abhyankar 
SNESVIProjectOntoBounds(SNES snes,Vec X)235d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X)
236d71ae5a4SJacob Faibussowitsch {
237c1a5e521SShri Abhyankar   const PetscScalar *xl, *xu;
238c1a5e521SShri Abhyankar   PetscScalar       *x;
239c1a5e521SShri Abhyankar   PetscInt           i, n;
240c1a5e521SShri Abhyankar 
241c1a5e521SShri Abhyankar   PetscFunctionBegin;
2429566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
2439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
2449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
2459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
246c1a5e521SShri Abhyankar 
247c1a5e521SShri Abhyankar   for (i = 0; i < n; i++) {
24810f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
24910f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
250c1a5e521SShri Abhyankar   }
2519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
2529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
2539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
2543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
255c1a5e521SShri Abhyankar }
256c1a5e521SShri Abhyankar 
257ceaaa498SBarry Smith /*@
258e4094ef1SJacob Faibussowitsch   SNESVIGetActiveSetIS - Gets the global indices for the active set variables
259693ddf92SShri Abhyankar 
260ceaaa498SBarry Smith   Input Parameters:
261ceaaa498SBarry Smith + snes - the `SNES` context
262ceaaa498SBarry Smith . X    - the `snes` solution vector
263e4094ef1SJacob Faibussowitsch - F    - the nonlinear function vector
264693ddf92SShri Abhyankar 
265e4094ef1SJacob Faibussowitsch   Output Parameter:
266693ddf92SShri Abhyankar . ISact - active set index set
267ceaaa498SBarry Smith 
268ceaaa498SBarry Smith   Level: developer
269ceaaa498SBarry Smith 
270420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`
271ceaaa498SBarry Smith @*/
SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS * ISact)272d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact)
273d71ae5a4SJacob Faibussowitsch {
274c2fc9fa9SBarry Smith   Vec                Xl = snes->xl, Xu = snes->xu;
275693ddf92SShri Abhyankar   const PetscScalar *x, *f, *xl, *xu;
276693ddf92SShri Abhyankar   PetscInt          *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0;
277349187a7SBarry Smith   PetscReal          zerotolerance = snes->vizerotolerance;
278d950fb63SShri Abhyankar 
279d950fb63SShri Abhyankar   PetscFunctionBegin;
2809566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nlocal));
2819566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh));
2829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xl, &xl));
2849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xu, &xu));
2859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
286693ddf92SShri Abhyankar   /* Compute active set size */
287d950fb63SShri Abhyankar   for (i = 0; i < nlocal; i++) {
288349187a7SBarry Smith     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++;
289d950fb63SShri Abhyankar   }
290693ddf92SShri Abhyankar 
2919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc_isact, &idx_act));
292d950fb63SShri Abhyankar 
293693ddf92SShri Abhyankar   /* Set active set indices */
294d950fb63SShri Abhyankar   for (i = 0; i < nlocal; i++) {
295349187a7SBarry Smith     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;
296d950fb63SShri Abhyankar   }
297d950fb63SShri Abhyankar 
298693ddf92SShri Abhyankar   /* Create active set IS */
2999566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact));
300d950fb63SShri Abhyankar 
3019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
3029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xl, &xl));
3039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xu, &xu));
3049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
306d950fb63SShri Abhyankar }
307d950fb63SShri Abhyankar 
308d5def619SJonas Heinzmann /*@
309d5def619SJonas Heinzmann   SNESVIComputeInactiveSetFnorm - Computes the function norm for variational inequalities on the inactive set
310d5def619SJonas Heinzmann 
311d5def619SJonas Heinzmann   Input Parameters:
312d5def619SJonas Heinzmann + snes - the `SNES` context
313d5def619SJonas Heinzmann . F    - the nonlinear function vector
314d5def619SJonas Heinzmann - X    - the `SNES` solution vector
315d5def619SJonas Heinzmann 
316d5def619SJonas Heinzmann   Output Parameter:
317d5def619SJonas Heinzmann . fnorm - the function norm
318d5def619SJonas Heinzmann 
319d5def619SJonas Heinzmann   Level: developer
320d5def619SJonas Heinzmann 
321d5def619SJonas Heinzmann .seealso: [](ch_snes), `SNES`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `SNESLineSearchSetVIFunctions()`
322d5def619SJonas Heinzmann @*/
SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X,PetscReal * fnorm)323d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm)
324d71ae5a4SJacob Faibussowitsch {
325fe835674SShri Abhyankar   const PetscScalar *x, *xl, *xu, *f;
326fe835674SShri Abhyankar   PetscInt           i, n;
327d5def619SJonas Heinzmann   PetscReal          zerotolerance = snes->vizerotolerance;
328fe835674SShri Abhyankar 
329fe835674SShri Abhyankar   PetscFunctionBegin;
330d5def619SJonas Heinzmann   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
331d5def619SJonas Heinzmann   PetscAssertPointer(fnorm, 4);
3329566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
3339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
3349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
3359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
3369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
337d5def619SJonas Heinzmann   *fnorm = 0.0;
338fe835674SShri Abhyankar   for (i = 0; i < n; i++) {
339d5def619SJonas Heinzmann     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]);
3408f918527SKarl Rupp   }
3419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
3429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
3439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
3449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
345d5def619SJonas Heinzmann   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
34662d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
348fe835674SShri Abhyankar }
349fe835674SShri Abhyankar 
350d5def619SJonas Heinzmann /*@
351d5def619SJonas Heinzmann   SNESVIComputeInactiveSetFtY - Computes the directional derivative for variational inequalities on the inactive set,
352d5def619SJonas Heinzmann   assuming that there exists some $G(x)$ for which the `SNESFunctionFn` $F(x) = grad G(x)$ (relevant for some line search algorithms)
353d5def619SJonas Heinzmann 
354d5def619SJonas Heinzmann   Input Parameters:
355d5def619SJonas Heinzmann + snes - the `SNES` context
356d5def619SJonas Heinzmann . F    - the nonlinear function vector
357d5def619SJonas Heinzmann . X    - the `SNES` solution vector
358d5def619SJonas Heinzmann - Y    - the direction vector
359d5def619SJonas Heinzmann 
360d5def619SJonas Heinzmann   Output Parameter:
361d5def619SJonas Heinzmann . fty - the directional derivative
362d5def619SJonas Heinzmann 
363d5def619SJonas Heinzmann   Level: developer
364d5def619SJonas Heinzmann 
365d5def619SJonas Heinzmann .seealso: [](ch_snes), `SNES`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`
366d5def619SJonas Heinzmann @*/
SNESVIComputeInactiveSetFtY(SNES snes,Vec F,Vec X,Vec Y,PetscScalar * fty)367d5def619SJonas Heinzmann PetscErrorCode SNESVIComputeInactiveSetFtY(SNES snes, Vec F, Vec X, Vec Y, PetscScalar *fty)
368d5def619SJonas Heinzmann {
369d5def619SJonas Heinzmann   const PetscScalar *x, *xl, *xu, *y, *f;
370d5def619SJonas Heinzmann   PetscInt           i, n;
371d5def619SJonas Heinzmann   PetscReal          zerotolerance = snes->vizerotolerance;
372d5def619SJonas Heinzmann 
373d5def619SJonas Heinzmann   PetscFunctionBegin;
374d5def619SJonas Heinzmann   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
375d5def619SJonas Heinzmann   PetscAssertPointer(fty, 5);
376d5def619SJonas Heinzmann   PetscCall(VecGetLocalSize(X, &n));
377d5def619SJonas Heinzmann   PetscCall(VecGetArrayRead(F, &f));
378d5def619SJonas Heinzmann   PetscCall(VecGetArrayRead(X, &x));
379d5def619SJonas Heinzmann   PetscCall(VecGetArrayRead(snes->xl, &xl));
380d5def619SJonas Heinzmann   PetscCall(VecGetArrayRead(snes->xu, &xu));
381d5def619SJonas Heinzmann   PetscCall(VecGetArrayRead(Y, &y));
382d5def619SJonas Heinzmann   *fty = 0.0;
383d5def619SJonas Heinzmann   for (i = 0; i < n; i++) {
384d5def619SJonas Heinzmann     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]);
385d5def619SJonas Heinzmann   }
386d5def619SJonas Heinzmann   PetscCall(VecRestoreArrayRead(F, &f));
387d5def619SJonas Heinzmann   PetscCall(VecRestoreArrayRead(X, &x));
388d5def619SJonas Heinzmann   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
389d5def619SJonas Heinzmann   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
390d5def619SJonas Heinzmann   PetscCall(VecRestoreArrayRead(Y, &y));
391d5def619SJonas Heinzmann   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, fty, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
392d5def619SJonas Heinzmann   PetscFunctionReturn(PETSC_SUCCESS);
393d5def619SJonas Heinzmann }
394d5def619SJonas Heinzmann 
SNESVIDMComputeVariableBounds(SNES snes,Vec xl,Vec xu)395ba38deedSJacob Faibussowitsch static PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu)
396d71ae5a4SJacob Faibussowitsch {
39708da532bSDmitry Karpeev   PetscFunctionBegin;
3989566063dSJacob Faibussowitsch   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
3993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40008da532bSDmitry Karpeev }
40108da532bSDmitry Karpeev 
4022b4ed282SShri Abhyankar /*
403c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
4042b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
4052b4ed282SShri Abhyankar 
4062b4ed282SShri Abhyankar    Input Parameter:
4072b4ed282SShri Abhyankar .  snes - the SNES context
4082b4ed282SShri Abhyankar 
4092b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
4102b4ed282SShri Abhyankar 
4112b4ed282SShri Abhyankar    Notes:
4122b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
4132b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
4142b4ed282SShri Abhyankar    the call to SNESSolve().
4152b4ed282SShri Abhyankar  */
SNESSetUp_VI(SNES snes)416d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_VI(SNES snes)
417d71ae5a4SJacob Faibussowitsch {
4182b4ed282SShri Abhyankar   PetscInt i_start[3], i_end[3];
4192b4ed282SShri Abhyankar 
4202b4ed282SShri Abhyankar   PetscFunctionBegin;
4219566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 1));
4229566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
4232b4ed282SShri Abhyankar 
42408da532bSDmitry Karpeev   if (!snes->ops->computevariablebounds && snes->dm) {
425a201590fSDmitry Karpeev     PetscBool flag;
4269566063dSJacob Faibussowitsch     PetscCall(DMHasVariableBounds(snes->dm, &flag));
427ad540459SPierre Jolivet     if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
42836b2a9f3SBarry Smith   }
429a201590fSDmitry Karpeev   if (!snes->usersetbounds) {
430c2fc9fa9SBarry Smith     if (snes->ops->computevariablebounds) {
431b2ccae6bSStefano Zampini       if (!snes->xl) PetscCall(VecDuplicate(snes->work[0], &snes->xl));
432b2ccae6bSStefano Zampini       if (!snes->xu) PetscCall(VecDuplicate(snes->work[0], &snes->xu));
433dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
4341aa26658SKarl Rupp     } else if (!snes->xl && !snes->xu) {
4352176524cSBarry Smith       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
436b2ccae6bSStefano Zampini       PetscCall(VecDuplicate(snes->work[0], &snes->xl));
4379566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xl, PETSC_NINFINITY));
438b2ccae6bSStefano Zampini       PetscCall(VecDuplicate(snes->work[0], &snes->xu));
4399566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xu, PETSC_INFINITY));
4402b4ed282SShri Abhyankar     } else {
4412b4ed282SShri Abhyankar       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
442b2ccae6bSStefano Zampini       PetscCall(VecGetOwnershipRange(snes->work[0], i_start, i_end));
4439566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1));
4449566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2));
4452b4ed282SShri Abhyankar       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]))
4462b4ed282SShri Abhyankar         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.");
4472b4ed282SShri Abhyankar     }
448a201590fSDmitry Karpeev   }
4493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4502b4ed282SShri Abhyankar }
SNESReset_VI(SNES snes)451d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_VI(SNES snes)
452d71ae5a4SJacob Faibussowitsch {
4532176524cSBarry Smith   PetscFunctionBegin;
4549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
4559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
4562d6615e8SDmitry Karpeev   snes->usersetbounds = PETSC_FALSE;
4573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4582176524cSBarry Smith }
4592176524cSBarry Smith 
4602b4ed282SShri Abhyankar /*
4612b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
4622b4ed282SShri Abhyankar    with SNESCreate_VI().
4632b4ed282SShri Abhyankar 
4642b4ed282SShri Abhyankar    Input Parameter:
4652b4ed282SShri Abhyankar .  snes - the SNES context
4662b4ed282SShri Abhyankar 
4672b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4682b4ed282SShri Abhyankar  */
SNESDestroy_VI(SNES snes)469d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_VI(SNES snes)
470d71ae5a4SJacob Faibussowitsch {
4712b4ed282SShri Abhyankar   PetscFunctionBegin;
4729566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
4732b4ed282SShri Abhyankar 
4742b4ed282SShri Abhyankar   /* clear composed functions */
4752e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL));
4762e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL));
4773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4782b4ed282SShri Abhyankar }
4797fe79bc4SShri Abhyankar 
4805559b345SBarry Smith /*@
481420bcc1bSBarry Smith   SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. `xl` <= x <= `xu`. This allows solving
482f6dfbefdSBarry Smith   (differential) variable inequalities.
4832b4ed282SShri Abhyankar 
4842b4ed282SShri Abhyankar   Input Parameters:
485f6dfbefdSBarry Smith + snes - the `SNES` context.
4862b4ed282SShri Abhyankar . xl   - lower bound.
487a2b725a8SWilliam Gropp - xu   - upper bound.
4882b4ed282SShri Abhyankar 
4892fe279fdSBarry Smith   Level: advanced
4902fe279fdSBarry Smith 
4912b4ed282SShri Abhyankar   Notes:
4922b4ed282SShri Abhyankar   If this routine is not called then the lower and upper bounds are set to
493f6dfbefdSBarry Smith   `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`.
49484c105d7SBarry Smith 
495420bcc1bSBarry Smith   Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS` or semi-smooth `SNESVINEWTONSSLS` solvers.
496f0b84518SBarry Smith 
497f6dfbefdSBarry Smith   For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY`
498f6dfbefdSBarry Smith 
49933e0caf2SBarry Smith   `SNESVISetComputeVariableBounds()` can be used to provide a function that computes the bounds. This should be used if you are using, for example, grid
500f6dfbefdSBarry Smith   sequencing and need bounds set for a variety of vectors
501f6dfbefdSBarry Smith 
5022913281cSPierre Jolivet .seealso: [](sec_vi), `SNES`, `SNESVIGetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
5035559b345SBarry Smith @*/
SNESVISetVariableBounds(SNES snes,Vec xl,Vec xu)504d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
505d71ae5a4SJacob Faibussowitsch {
5065f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES, Vec, Vec);
5072b4ed282SShri Abhyankar 
5082b4ed282SShri Abhyankar   PetscFunctionBegin;
5092b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
5102b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl, VEC_CLASSID, 2);
5112b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu, VEC_CLASSID, 3);
5129566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f));
513cac4c232SBarry Smith   if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu));
5149566063dSJacob Faibussowitsch   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
515a201590fSDmitry Karpeev   snes->usersetbounds = PETSC_TRUE;
5163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51761589011SJed Brown }
51861589011SJed Brown 
SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)519d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu)
520d71ae5a4SJacob Faibussowitsch {
52161589011SJed Brown   const PetscScalar *xxl, *xxu;
52261589011SJed Brown   PetscInt           i, n, cnt = 0;
52361589011SJed Brown 
52461589011SJed Brown   PetscFunctionBegin;
5259566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL));
5265f80ce2aSJacob Faibussowitsch   PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first");
527077aedafSJed Brown   {
528077aedafSJed Brown     PetscInt xlN, xuN, N;
5299566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xl, &xlN));
5309566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xu, &xuN));
5319566063dSJacob Faibussowitsch     PetscCall(VecGetSize(snes->vec_func, &N));
53263a3b9bcSJacob Faibussowitsch     PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N);
53363a3b9bcSJacob Faibussowitsch     PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N);
534077aedafSJed Brown   }
5359566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xl));
5369566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xu));
5379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
5389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
539c2fc9fa9SBarry Smith   snes->xl = xl;
540c2fc9fa9SBarry Smith   snes->xu = xu;
5419566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xl, &n));
5429566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xl, &xxl));
5439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xu, &xxu));
544e270355aSBarry Smith   for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
5451aa26658SKarl Rupp 
546462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
5479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xl, &xxl));
5489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xu, &xxu));
5493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5502b4ed282SShri Abhyankar }
55192c02d66SPeter Brune 
552cf836535SBlaise Bourdin /*@
553420bcc1bSBarry Smith   SNESVIGetVariableBounds - Gets the lower and upper bounds for the solution vector. `xl` <= x <= `xu`. These are used in solving
554cf836535SBlaise Bourdin   (differential) variable inequalities.
555cf836535SBlaise Bourdin 
556cf836535SBlaise Bourdin   Input Parameters:
557cf836535SBlaise Bourdin + snes - the `SNES` context.
558cf836535SBlaise Bourdin . xl   - lower bound (may be `NULL`)
559cf836535SBlaise Bourdin - xu   - upper bound (may be `NULL`)
560cf836535SBlaise Bourdin 
5612fe279fdSBarry Smith   Level: advanced
5622fe279fdSBarry Smith 
563420bcc1bSBarry Smith   Note:
564cf836535SBlaise Bourdin   These vectors are owned by the `SNESVI` and should not be destroyed by the caller
565cf836535SBlaise Bourdin 
5662913281cSPierre Jolivet .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
567cf836535SBlaise Bourdin @*/
SNESVIGetVariableBounds(SNES snes,Vec * xl,Vec * xu)568cf836535SBlaise Bourdin PetscErrorCode SNESVIGetVariableBounds(SNES snes, Vec *xl, Vec *xu)
569cf836535SBlaise Bourdin {
570cf836535SBlaise Bourdin   PetscFunctionBegin;
571cf836535SBlaise Bourdin   PetscCheck(snes->usersetbounds, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must set SNESVI bounds before calling SNESVIGetVariableBounds()");
572cf836535SBlaise Bourdin   if (xl) *xl = snes->xl;
573cf836535SBlaise Bourdin   if (xu) *xu = snes->xu;
5743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
575cf836535SBlaise Bourdin }
576cf836535SBlaise Bourdin 
SNESSetFromOptions_VI(SNES snes,PetscOptionItems PetscOptionsObject)577ce78bad3SBarry Smith PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems PetscOptionsObject)
578d71ae5a4SJacob Faibussowitsch {
5798afaa268SBarry Smith   PetscBool flg = PETSC_FALSE;
5802b4ed282SShri Abhyankar 
5812b4ed282SShri Abhyankar   PetscFunctionBegin;
582d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options");
5839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL));
5849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL));
5851baa6e33SBarry Smith   if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL));
586de34d3e9SBarry Smith   flg = PETSC_FALSE;
587b6f515dbSMatthew G. Knepley   PetscCall(SNESMonitorSetFromOptions(snes, "-snes_vi_monitor_residual", "View residual at each iteration, using zero for active constraints", "SNESVIMonitorResidual", SNESVIMonitorResidual, NULL));
588b6f515dbSMatthew G. Knepley   PetscCall(SNESMonitorSetFromOptions(snes, "-snes_vi_monitor_active", "View active set at each iteration, using zero for inactive dofs", "SNESVIMonitorActive", SNESVIMonitorActive, NULL));
589d0609cedSBarry Smith   PetscOptionsHeadEnd();
5903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
591ed70e96aSJungho Lee }
592