xref: /petsc/src/snes/impls/vi/vi.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
21e25c274SJed Brown #include <petscdm.h>
3d0af7cd3SBarry Smith 
42176524cSBarry Smith /*@C
5*f6dfbefdSBarry Smith    SNESVISetComputeVariableBounds - Sets a function that is called to compute the bounds on variable for
6*f6dfbefdSBarry Smith    (differential) variable inequalities.
72176524cSBarry Smith 
87a7aea1fSJed Brown    Input parameter:
9*f6dfbefdSBarry Smith +  snes - the `SNES` context
10*f6dfbefdSBarry Smith -  compute - function that computes the bounds
11*f6dfbefdSBarry Smith 
12*f6dfbefdSBarry Smith  Calling Sequence of function:
13*f6dfbefdSBarry Smith   PetscErrorCode compute(SNES snes,Vec lower,Vec higher, void *ctx)
14*f6dfbefdSBarry Smith 
15*f6dfbefdSBarry Smith + snes - the `SNES` context
16*f6dfbefdSBarry Smith . lower - vector to hold lower bounds
17*f6dfbefdSBarry Smith - higher - vector to hold upper bounds
182176524cSBarry Smith 
192bd2b0e6SSatish Balay    Level: advanced
202176524cSBarry Smith 
21*f6dfbefdSBarry Smith    Notes:
22f0b84518SBarry Smith    Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.
23f0b84518SBarry Smith 
24*f6dfbefdSBarry Smith    For entries with no bounds you can set `PETSC_NINFINITY` or `PETSC_INFINITY`
25c1c3a0ecSBarry Smith 
26*f6dfbefdSBarry Smith    You may use `SNESVISetVariableBounds()` to provide the bounds once if they will never change
27*f6dfbefdSBarry Smith 
28*f6dfbefdSBarry Smith    If you have associated a `DM` with the `SNES` and provided a function to the `DM` via `DMSetVariableBounds()` that will be used automatically
29*f6dfbefdSBarry Smith    to provide the bounds and you need not use this function.
30*f6dfbefdSBarry Smith 
31*f6dfbefdSBarry Smith .seealso: `SNESVISetVariableBounds()`, `DMSetVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`,
32*f6dfbefdSBarry Smith           'SNESSetType()`
33aab9d709SJed Brown @*/
349371c9d4SSatish Balay PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES, Vec, Vec)) {
355f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES, PetscErrorCode(*)(SNES, Vec, Vec));
362176524cSBarry Smith 
372176524cSBarry Smith   PetscFunctionBegin;
3861589011SJed Brown   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
399566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", &f));
40cac4c232SBarry Smith   if (f) PetscUseMethod(snes, "SNESVISetComputeVariableBounds_C", (SNES, PetscErrorCode(*)(SNES, Vec, Vec)), (snes, compute));
419566063dSJacob Faibussowitsch   else PetscCall(SNESVISetComputeVariableBounds_VI(snes, compute));
422176524cSBarry Smith   PetscFunctionReturn(0);
432176524cSBarry Smith }
442176524cSBarry Smith 
459371c9d4SSatish Balay PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes, SNESVIComputeVariableBoundsFunction compute) {
4661589011SJed Brown   PetscFunctionBegin;
4761589011SJed Brown   snes->ops->computevariablebounds = compute;
4861589011SJed Brown   PetscFunctionReturn(0);
4961589011SJed Brown }
502176524cSBarry Smith 
519371c9d4SSatish Balay PetscErrorCode SNESVIMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy) {
52ffdf2a20SBarry Smith   Vec         X, F, Finactive;
53ffdf2a20SBarry Smith   IS          isactive;
54ffdf2a20SBarry Smith   PetscViewer viewer = (PetscViewer)dummy;
55ffdf2a20SBarry Smith 
56ffdf2a20SBarry Smith   PetscFunctionBegin;
579566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
589566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &X));
599566063dSJacob Faibussowitsch   PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive));
609566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(F, &Finactive));
619566063dSJacob Faibussowitsch   PetscCall(VecCopy(F, Finactive));
629566063dSJacob Faibussowitsch   PetscCall(VecISSet(Finactive, isactive, 0.0));
639566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isactive));
649566063dSJacob Faibussowitsch   PetscCall(VecView(Finactive, viewer));
659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Finactive));
66ffdf2a20SBarry Smith   PetscFunctionReturn(0);
67ffdf2a20SBarry Smith }
68ffdf2a20SBarry Smith 
699371c9d4SSatish Balay PetscErrorCode SNESMonitorVI(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy) {
704d4332d5SBarry Smith   PetscViewer        viewer = (PetscViewer)dummy;
719308c137SBarry Smith   const PetscScalar *x, *xl, *xu, *f;
726fd67740SBarry Smith   PetscInt           i, n, act[2] = {0, 0}, fact[2], N;
736a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
746a9e2d46SJungho Lee   PetscInt           act_bound[2] = {0, 0}, fact_bound[2];
75349187a7SBarry Smith   PetscReal          rnorm, fnorm, zerotolerance = snes->vizerotolerance;
769d1809e2SSatish Balay   double             tmp;
779308c137SBarry Smith 
789308c137SBarry Smith   PetscFunctionBegin;
794d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
809566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
819566063dSJacob Faibussowitsch   PetscCall(VecGetSize(snes->vec_sol, &N));
829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->vec_sol, &x));
859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->vec_func, &f));
869308c137SBarry Smith 
879308c137SBarry Smith   rnorm = 0.0;
889308c137SBarry Smith   for (i = 0; i < n; i++) {
89349187a7SBarry 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))) rnorm += PetscRealPart(PetscConj(f[i]) * f[i]);
90349187a7SBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++;
91349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++;
92ce94432eSBarry Smith     else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Can never get here");
939308c137SBarry Smith   }
946a9e2d46SJungho Lee 
956a9e2d46SJungho Lee   for (i = 0; i < n; i++) {
96349187a7SBarry Smith     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++;
97349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++;
986a9e2d46SJungho Lee   }
999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->vec_func, &f));
1009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
1019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
1029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->vec_sol, &x));
1031c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&rnorm, &fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
1041c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(act, fact, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
1051c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(act_bound, fact_bound, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
106f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1076fd67740SBarry Smith 
1089566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
1099d1809e2SSatish Balay   if (snes->ntruebounds) tmp = ((double)(fact[0] + fact[1])) / ((double)snes->ntruebounds);
1109d1809e2SSatish Balay   else tmp = 0.0;
11163a3b9bcSJacob 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));
1126a9e2d46SJungho Lee 
1139566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
1149308c137SBarry Smith   PetscFunctionReturn(0);
1159308c137SBarry Smith }
1169308c137SBarry Smith 
1172b4ed282SShri Abhyankar /*
1182b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1192b4ed282SShri 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
1202b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1212b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1222b4ed282SShri Abhyankar */
1239371c9d4SSatish Balay PetscErrorCode SNESVICheckLocalMin_Private(SNES snes, Mat A, Vec F, Vec W, PetscReal fnorm, PetscBool *ismin) {
1242b4ed282SShri Abhyankar   PetscReal a1;
125ace3abfcSBarry Smith   PetscBool hastranspose;
1262b4ed282SShri Abhyankar 
1272b4ed282SShri Abhyankar   PetscFunctionBegin;
1282b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1299566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
1302b4ed282SShri Abhyankar   if (hastranspose) {
1312b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1329566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, F, W));
1339566063dSJacob Faibussowitsch     PetscCall(VecNorm(W, NORM_2, &a1));
1349566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "|| J^T F|| %g near zero implies found a local minimum\n", (double)(a1 / fnorm)));
1352b4ed282SShri Abhyankar     if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE;
1362b4ed282SShri Abhyankar   } else {
1372b4ed282SShri Abhyankar     Vec         work;
1382b4ed282SShri Abhyankar     PetscScalar result;
1392b4ed282SShri Abhyankar     PetscReal   wnorm;
1402b4ed282SShri Abhyankar 
1419566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(W, NULL));
1429566063dSJacob Faibussowitsch     PetscCall(VecNorm(W, NORM_2, &wnorm));
1439566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(W, &work));
1449566063dSJacob Faibussowitsch     PetscCall(MatMult(A, W, work));
1459566063dSJacob Faibussowitsch     PetscCall(VecDot(F, work, &result));
1469566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&work));
1472b4ed282SShri Abhyankar     a1 = PetscAbsScalar(result) / (fnorm * wnorm);
1489566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n", (double)a1));
1492b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1502b4ed282SShri Abhyankar   }
1512b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1522b4ed282SShri Abhyankar }
1532b4ed282SShri Abhyankar 
1542b4ed282SShri Abhyankar /*
1552b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
1562b4ed282SShri Abhyankar */
1579371c9d4SSatish Balay PetscErrorCode SNESVICheckResidual_Private(SNES snes, Mat A, Vec F, Vec X, Vec W1, Vec W2) {
1582b4ed282SShri Abhyankar   PetscReal a1, a2;
159ace3abfcSBarry Smith   PetscBool hastranspose;
1602b4ed282SShri Abhyankar 
1612b4ed282SShri Abhyankar   PetscFunctionBegin;
1629566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
1632b4ed282SShri Abhyankar   if (hastranspose) {
1649566063dSJacob Faibussowitsch     PetscCall(MatMult(A, X, W1));
1659566063dSJacob Faibussowitsch     PetscCall(VecAXPY(W1, -1.0, F));
1662b4ed282SShri Abhyankar 
1672b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1689566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, W1, W2));
1699566063dSJacob Faibussowitsch     PetscCall(VecNorm(W1, NORM_2, &a1));
1709566063dSJacob Faibussowitsch     PetscCall(VecNorm(W2, NORM_2, &a2));
17148a46eb9SPierre Jolivet     if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n", (double)(a2 / a1)));
1722b4ed282SShri Abhyankar   }
1732b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1742b4ed282SShri Abhyankar }
1752b4ed282SShri Abhyankar 
1762b4ed282SShri Abhyankar /*
1778d359177SBarry Smith   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
1782b4ed282SShri Abhyankar 
1792b4ed282SShri Abhyankar   Notes:
1802b4ed282SShri Abhyankar   The convergence criterion currently implemented is
1812b4ed282SShri Abhyankar   merit < abstol
1822b4ed282SShri Abhyankar   merit < rtol*merit_initial
1832b4ed282SShri Abhyankar */
1849371c9d4SSatish Balay PetscErrorCode SNESConvergedDefault_VI(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gradnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) {
1852b4ed282SShri Abhyankar   PetscFunctionBegin;
1862b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1872b4ed282SShri Abhyankar   PetscValidPointer(reason, 6);
1882b4ed282SShri Abhyankar 
1892b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
1902b4ed282SShri Abhyankar 
1912b4ed282SShri Abhyankar   if (!it) {
1922b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
1937fe79bc4SShri Abhyankar     snes->ttol = fnorm * snes->rtol;
1942b4ed282SShri Abhyankar   }
1957fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
1969566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n"));
1972b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
1980d6f27a8SBarry Smith   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
1999566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g\n", (double)fnorm, (double)snes->abstol));
2002b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
201e71169deSBarry Smith   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
20263a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs));
2032b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
2042b4ed282SShri Abhyankar   }
2052b4ed282SShri Abhyankar 
2062b4ed282SShri Abhyankar   if (it && !*reason) {
2077fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2089566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g (relative tolerance)\n", (double)fnorm, (double)snes->ttol));
2092b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2102b4ed282SShri Abhyankar     }
2112b4ed282SShri Abhyankar   }
2122b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2132b4ed282SShri Abhyankar }
2142b4ed282SShri Abhyankar 
215c1a5e521SShri Abhyankar /*
216c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
217c1a5e521SShri Abhyankar 
218c1a5e521SShri Abhyankar    Input Parameters:
219c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
220c1a5e521SShri Abhyankar 
221c1a5e521SShri Abhyankar    Output Parameters:
222c1a5e521SShri Abhyankar .  X - Bound projected X
223c1a5e521SShri Abhyankar 
224c1a5e521SShri Abhyankar */
225c1a5e521SShri Abhyankar 
2269371c9d4SSatish Balay PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X) {
227c1a5e521SShri Abhyankar   const PetscScalar *xl, *xu;
228c1a5e521SShri Abhyankar   PetscScalar       *x;
229c1a5e521SShri Abhyankar   PetscInt           i, n;
230c1a5e521SShri Abhyankar 
231c1a5e521SShri Abhyankar   PetscFunctionBegin;
2329566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
2339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
2349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
2359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
236c1a5e521SShri Abhyankar 
237c1a5e521SShri Abhyankar   for (i = 0; i < n; i++) {
23810f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
23910f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
240c1a5e521SShri Abhyankar   }
2419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
2429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
2439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
244c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
245c1a5e521SShri Abhyankar }
246c1a5e521SShri Abhyankar 
247693ddf92SShri Abhyankar /*
248693ddf92SShri Abhyankar    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
249693ddf92SShri Abhyankar 
2507a7aea1fSJed Brown    Input parameter:
251693ddf92SShri Abhyankar .  snes - the SNES context
252693ddf92SShri Abhyankar .  X    - the snes solution vector
253693ddf92SShri Abhyankar .  F    - the nonlinear function vector
254693ddf92SShri Abhyankar 
2557a7aea1fSJed Brown    Output parameter:
256693ddf92SShri Abhyankar .  ISact - active set index set
257693ddf92SShri Abhyankar  */
2589371c9d4SSatish Balay PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact) {
259c2fc9fa9SBarry Smith   Vec                Xl = snes->xl, Xu = snes->xu;
260693ddf92SShri Abhyankar   const PetscScalar *x, *f, *xl, *xu;
261693ddf92SShri Abhyankar   PetscInt          *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0;
262349187a7SBarry Smith   PetscReal          zerotolerance = snes->vizerotolerance;
263d950fb63SShri Abhyankar 
264d950fb63SShri Abhyankar   PetscFunctionBegin;
2659566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nlocal));
2669566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh));
2679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xl, &xl));
2699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xu, &xu));
2709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
271693ddf92SShri Abhyankar   /* Compute active set size */
272d950fb63SShri Abhyankar   for (i = 0; i < nlocal; i++) {
273349187a7SBarry 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++;
274d950fb63SShri Abhyankar   }
275693ddf92SShri Abhyankar 
2769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc_isact, &idx_act));
277d950fb63SShri Abhyankar 
278693ddf92SShri Abhyankar   /* Set active set indices */
279d950fb63SShri Abhyankar   for (i = 0; i < nlocal; i++) {
280349187a7SBarry 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;
281d950fb63SShri Abhyankar   }
282d950fb63SShri Abhyankar 
283693ddf92SShri Abhyankar   /* Create active set IS */
2849566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact));
285d950fb63SShri Abhyankar 
2869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xl, &xl));
2889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xu, &xu));
2899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
290d950fb63SShri Abhyankar   PetscFunctionReturn(0);
291d950fb63SShri Abhyankar }
292d950fb63SShri Abhyankar 
2939371c9d4SSatish Balay PetscErrorCode SNESVICreateIndexSets_RS(SNES snes, Vec X, Vec F, IS *ISact, IS *ISinact) {
294077aedafSJed Brown   PetscInt rstart, rend;
295693ddf92SShri Abhyankar 
296693ddf92SShri Abhyankar   PetscFunctionBegin;
2979566063dSJacob Faibussowitsch   PetscCall(SNESVIGetActiveSetIS(snes, X, F, ISact));
2989566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &rstart, &rend));
2999566063dSJacob Faibussowitsch   PetscCall(ISComplement(*ISact, rstart, rend, ISinact));
300693ddf92SShri Abhyankar   PetscFunctionReturn(0);
301693ddf92SShri Abhyankar }
302693ddf92SShri Abhyankar 
3039371c9d4SSatish Balay PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm) {
304fe835674SShri Abhyankar   const PetscScalar *x, *xl, *xu, *f;
305fe835674SShri Abhyankar   PetscInt           i, n;
306feb237baSPierre Jolivet   PetscReal          rnorm, zerotolerance = snes->vizerotolerance;
307fe835674SShri Abhyankar 
308fe835674SShri Abhyankar   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
3109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
3119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
3129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
3139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
314fe835674SShri Abhyankar   rnorm = 0.0;
315fe835674SShri Abhyankar   for (i = 0; i < n; i++) {
316349187a7SBarry 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))) rnorm += PetscRealPart(PetscConj(f[i]) * f[i]);
3178f918527SKarl Rupp   }
3189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
3199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
3209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
3219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
3221c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
32362d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
324fe835674SShri Abhyankar   PetscFunctionReturn(0);
325fe835674SShri Abhyankar }
326fe835674SShri Abhyankar 
3279371c9d4SSatish Balay PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu) {
32808da532bSDmitry Karpeev   PetscFunctionBegin;
3299566063dSJacob Faibussowitsch   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
33008da532bSDmitry Karpeev   PetscFunctionReturn(0);
33108da532bSDmitry Karpeev }
33208da532bSDmitry Karpeev 
3332b4ed282SShri Abhyankar /*
334c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3352b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3362b4ed282SShri Abhyankar 
3372b4ed282SShri Abhyankar    Input Parameter:
3382b4ed282SShri Abhyankar .  snes - the SNES context
3392b4ed282SShri Abhyankar 
3402b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
3412b4ed282SShri Abhyankar 
3422b4ed282SShri Abhyankar    Notes:
3432b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
3442b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
3452b4ed282SShri Abhyankar    the call to SNESSolve().
3462b4ed282SShri Abhyankar  */
3479371c9d4SSatish Balay PetscErrorCode SNESSetUp_VI(SNES snes) {
3482b4ed282SShri Abhyankar   PetscInt i_start[3], i_end[3];
3492b4ed282SShri Abhyankar 
3502b4ed282SShri Abhyankar   PetscFunctionBegin;
3519566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 1));
3529566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
3532b4ed282SShri Abhyankar 
35408da532bSDmitry Karpeev   if (!snes->ops->computevariablebounds && snes->dm) {
355a201590fSDmitry Karpeev     PetscBool flag;
3569566063dSJacob Faibussowitsch     PetscCall(DMHasVariableBounds(snes->dm, &flag));
357ad540459SPierre Jolivet     if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
35836b2a9f3SBarry Smith   }
359a201590fSDmitry Karpeev   if (!snes->usersetbounds) {
360c2fc9fa9SBarry Smith     if (snes->ops->computevariablebounds) {
3619566063dSJacob Faibussowitsch       if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3629566063dSJacob Faibussowitsch       if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
363dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
3641aa26658SKarl Rupp     } else if (!snes->xl && !snes->xu) {
3652176524cSBarry Smith       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
3669566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3679566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xl, PETSC_NINFINITY));
3689566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
3699566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xu, PETSC_INFINITY));
3702b4ed282SShri Abhyankar     } else {
3712b4ed282SShri Abhyankar       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
3729566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end));
3739566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1));
3749566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2));
3752b4ed282SShri 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]))
3762b4ed282SShri 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.");
3772b4ed282SShri Abhyankar     }
378a201590fSDmitry Karpeev   }
3792b4ed282SShri Abhyankar   PetscFunctionReturn(0);
3802b4ed282SShri Abhyankar }
3819371c9d4SSatish Balay PetscErrorCode SNESReset_VI(SNES snes) {
3822176524cSBarry Smith   PetscFunctionBegin;
3839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
3849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
3852d6615e8SDmitry Karpeev   snes->usersetbounds = PETSC_FALSE;
3862176524cSBarry Smith   PetscFunctionReturn(0);
3872176524cSBarry Smith }
3882176524cSBarry Smith 
3892b4ed282SShri Abhyankar /*
3902b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
3912b4ed282SShri Abhyankar    with SNESCreate_VI().
3922b4ed282SShri Abhyankar 
3932b4ed282SShri Abhyankar    Input Parameter:
3942b4ed282SShri Abhyankar .  snes - the SNES context
3952b4ed282SShri Abhyankar 
3962b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
3972b4ed282SShri Abhyankar  */
3989371c9d4SSatish Balay PetscErrorCode SNESDestroy_VI(SNES snes) {
3992b4ed282SShri Abhyankar   PetscFunctionBegin;
4009566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
4012b4ed282SShri Abhyankar 
4022b4ed282SShri Abhyankar   /* clear composed functions */
4032e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL));
4042e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL));
4052b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4062b4ed282SShri Abhyankar }
4077fe79bc4SShri Abhyankar 
4085559b345SBarry Smith /*@
409f0b84518SBarry Smith    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving
410*f6dfbefdSBarry Smith    (differential) variable inequalities.
4112b4ed282SShri Abhyankar 
4122b4ed282SShri Abhyankar    Input Parameters:
413*f6dfbefdSBarry Smith +  snes - the `SNES` context.
4142b4ed282SShri Abhyankar .  xl   - lower bound.
415a2b725a8SWilliam Gropp -  xu   - upper bound.
4162b4ed282SShri Abhyankar 
4172b4ed282SShri Abhyankar    Notes:
4182b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
419*f6dfbefdSBarry Smith    `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`.
42084c105d7SBarry Smith 
421f0b84518SBarry Smith    Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.
422f0b84518SBarry Smith 
423*f6dfbefdSBarry Smith    For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY`
424*f6dfbefdSBarry Smith 
425*f6dfbefdSBarry Smith    `SNESVISetVariableBounds()` can be used to provide a function that computes the bounds. This should be used if you are using, for example, grid
426*f6dfbefdSBarry Smith    sequencing and need bounds set for a variety of vectors
427*f6dfbefdSBarry Smith 
4282bd2b0e6SSatish Balay    Level: advanced
4292bd2b0e6SSatish Balay 
430f0b84518SBarry Smith .seealso: `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, SNESVINEWTONRSLS, SNESVINEWTONSSLS, 'SNESSetType()`
4315559b345SBarry Smith @*/
4329371c9d4SSatish Balay PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) {
4335f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES, Vec, Vec);
4342b4ed282SShri Abhyankar 
4352b4ed282SShri Abhyankar   PetscFunctionBegin;
4362b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
4372b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl, VEC_CLASSID, 2);
4382b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu, VEC_CLASSID, 3);
4399566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f));
440cac4c232SBarry Smith   if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu));
4419566063dSJacob Faibussowitsch   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
442a201590fSDmitry Karpeev   snes->usersetbounds = PETSC_TRUE;
44361589011SJed Brown   PetscFunctionReturn(0);
44461589011SJed Brown }
44561589011SJed Brown 
4469371c9d4SSatish Balay PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu) {
44761589011SJed Brown   const PetscScalar *xxl, *xxu;
44861589011SJed Brown   PetscInt           i, n, cnt = 0;
44961589011SJed Brown 
45061589011SJed Brown   PetscFunctionBegin;
4519566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL));
4525f80ce2aSJacob Faibussowitsch   PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first");
453077aedafSJed Brown   {
454077aedafSJed Brown     PetscInt xlN, xuN, N;
4559566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xl, &xlN));
4569566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xu, &xuN));
4579566063dSJacob Faibussowitsch     PetscCall(VecGetSize(snes->vec_func, &N));
45863a3b9bcSJacob Faibussowitsch     PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N);
45963a3b9bcSJacob Faibussowitsch     PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N);
460077aedafSJed Brown   }
4619566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xl));
4629566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xu));
4639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
4649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
465c2fc9fa9SBarry Smith   snes->xl = xl;
466c2fc9fa9SBarry Smith   snes->xu = xu;
4679566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xl, &n));
4689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xl, &xxl));
4699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xu, &xxu));
470e270355aSBarry Smith   for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
4711aa26658SKarl Rupp 
4721c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
4739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xl, &xxl));
4749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xu, &xxu));
4752b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4762b4ed282SShri Abhyankar }
47792c02d66SPeter Brune 
4789371c9d4SSatish Balay PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject) {
4798afaa268SBarry Smith   PetscBool flg = PETSC_FALSE;
4802b4ed282SShri Abhyankar 
4812b4ed282SShri Abhyankar   PetscFunctionBegin;
482d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options");
4839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL));
4849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL));
4851baa6e33SBarry Smith   if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL));
486de34d3e9SBarry Smith   flg = PETSC_FALSE;
4879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor_residual", "Monitor residual all non-active variables; using zero for active constraints", "SNESMonitorVIResidual", flg, &flg, NULL));
4881baa6e33SBarry Smith   if (flg) PetscCall(SNESMonitorSet(snes, SNESVIMonitorResidual, PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)), NULL));
489d0609cedSBarry Smith   PetscOptionsHeadEnd();
490ed70e96aSJungho Lee   PetscFunctionReturn(0);
491ed70e96aSJungho Lee }
492