xref: /petsc/src/snes/impls/vi/vi.c (revision 945b2ebdc5b7ebcb622acfbe60dbc4e0db537728)
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`,
30*945b2ebdSBarry Smith           `'SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
31aab9d709SJed Brown @*/
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 
44d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes, SNESVIComputeVariableBoundsFunction compute)
45d71ae5a4SJacob Faibussowitsch {
4661589011SJed Brown   PetscFunctionBegin;
4761589011SJed Brown   snes->ops->computevariablebounds = compute;
483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4961589011SJed Brown }
502176524cSBarry Smith 
51ba38deedSJacob Faibussowitsch static PetscErrorCode SNESVIMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy)
52d71ae5a4SJacob Faibussowitsch {
53ffdf2a20SBarry Smith   Vec         X, F, Finactive;
54ffdf2a20SBarry Smith   IS          isactive;
55ffdf2a20SBarry Smith   PetscViewer viewer = (PetscViewer)dummy;
56ffdf2a20SBarry Smith 
57ffdf2a20SBarry Smith   PetscFunctionBegin;
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));
629566063dSJacob Faibussowitsch   PetscCall(VecCopy(F, Finactive));
639566063dSJacob Faibussowitsch   PetscCall(VecISSet(Finactive, isactive, 0.0));
649566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isactive));
659566063dSJacob Faibussowitsch   PetscCall(VecView(Finactive, viewer));
669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Finactive));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68ffdf2a20SBarry Smith }
69ffdf2a20SBarry Smith 
70ba38deedSJacob Faibussowitsch static PetscErrorCode SNESMonitorVI(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy)
71d71ae5a4SJacob Faibussowitsch {
724d4332d5SBarry Smith   PetscViewer        viewer = (PetscViewer)dummy;
739308c137SBarry Smith   const PetscScalar *x, *xl, *xu, *f;
746fd67740SBarry Smith   PetscInt           i, n, act[2] = {0, 0}, fact[2], N;
756a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
766a9e2d46SJungho Lee   PetscInt  act_bound[2] = {0, 0}, fact_bound[2];
77349187a7SBarry Smith   PetscReal rnorm, fnorm, zerotolerance = snes->vizerotolerance;
789d1809e2SSatish Balay   double    tmp;
799308c137SBarry Smith 
809308c137SBarry Smith   PetscFunctionBegin;
814d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
829566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
839566063dSJacob Faibussowitsch   PetscCall(VecGetSize(snes->vec_sol, &N));
849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->vec_sol, &x));
879566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->vec_func, &f));
889308c137SBarry Smith 
899308c137SBarry Smith   rnorm = 0.0;
909308c137SBarry Smith   for (i = 0; i < n; i++) {
91349187a7SBarry 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]);
92349187a7SBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++;
93349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++;
94ce94432eSBarry Smith     else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Can never get here");
959308c137SBarry Smith   }
966a9e2d46SJungho Lee 
976a9e2d46SJungho Lee   for (i = 0; i < n; i++) {
98349187a7SBarry Smith     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++;
99349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++;
1006a9e2d46SJungho Lee   }
1019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->vec_func, &f));
1029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
1039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
1049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->vec_sol, &x));
1051c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&rnorm, &fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
1061c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(act, fact, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
1071c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(act_bound, fact_bound, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
108f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1096fd67740SBarry Smith 
1109566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
1119d1809e2SSatish Balay   if (snes->ntruebounds) tmp = ((double)(fact[0] + fact[1])) / ((double)snes->ntruebounds);
1129d1809e2SSatish Balay   else tmp = 0.0;
11363a3b9bcSJacob 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));
1146a9e2d46SJungho Lee 
1159566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
1163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1179308c137SBarry Smith }
1189308c137SBarry Smith 
1192b4ed282SShri Abhyankar /*
1202b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1212b4ed282SShri 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
1222b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1232b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1242b4ed282SShri Abhyankar */
125d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVICheckLocalMin_Private(SNES snes, Mat A, Vec F, Vec W, PetscReal fnorm, PetscBool *ismin)
126d71ae5a4SJacob Faibussowitsch {
1272b4ed282SShri Abhyankar   PetscReal a1;
128ace3abfcSBarry Smith   PetscBool hastranspose;
1292b4ed282SShri Abhyankar 
1302b4ed282SShri Abhyankar   PetscFunctionBegin;
1312b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1329566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
1332b4ed282SShri Abhyankar   if (hastranspose) {
1342b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1359566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, F, W));
1369566063dSJacob Faibussowitsch     PetscCall(VecNorm(W, NORM_2, &a1));
1379566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "|| J^T F|| %g near zero implies found a local minimum\n", (double)(a1 / fnorm)));
1382b4ed282SShri Abhyankar     if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE;
1392b4ed282SShri Abhyankar   } else {
1402b4ed282SShri Abhyankar     Vec         work;
1412b4ed282SShri Abhyankar     PetscScalar result;
1422b4ed282SShri Abhyankar     PetscReal   wnorm;
1432b4ed282SShri Abhyankar 
1449566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(W, NULL));
1459566063dSJacob Faibussowitsch     PetscCall(VecNorm(W, NORM_2, &wnorm));
1469566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(W, &work));
1479566063dSJacob Faibussowitsch     PetscCall(MatMult(A, W, work));
1489566063dSJacob Faibussowitsch     PetscCall(VecDot(F, work, &result));
1499566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&work));
1502b4ed282SShri Abhyankar     a1 = PetscAbsScalar(result) / (fnorm * wnorm);
1519566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n", (double)a1));
1522b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1532b4ed282SShri Abhyankar   }
1543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1552b4ed282SShri Abhyankar }
1562b4ed282SShri Abhyankar 
1572b4ed282SShri Abhyankar /*
1588d359177SBarry Smith   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
1592b4ed282SShri Abhyankar 
1602b4ed282SShri Abhyankar   Notes:
1612b4ed282SShri Abhyankar   The convergence criterion currently implemented is
1622b4ed282SShri Abhyankar   merit < abstol
1632b4ed282SShri Abhyankar   merit < rtol*merit_initial
1642b4ed282SShri Abhyankar */
165d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESConvergedDefault_VI(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gradnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
166d71ae5a4SJacob Faibussowitsch {
1672b4ed282SShri Abhyankar   PetscFunctionBegin;
1682b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1694f572ea9SToby Isaac   PetscAssertPointer(reason, 6);
1702b4ed282SShri Abhyankar 
1712b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
1722b4ed282SShri Abhyankar 
1732b4ed282SShri Abhyankar   if (!it) {
1742b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
1757fe79bc4SShri Abhyankar     snes->ttol = fnorm * snes->rtol;
1762b4ed282SShri Abhyankar   }
1777fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
1789566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n"));
1792b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
1800d6f27a8SBarry Smith   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
1819566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g\n", (double)fnorm, (double)snes->abstol));
1822b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
183e71169deSBarry Smith   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
18463a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs));
1852b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
1862b4ed282SShri Abhyankar   }
1872b4ed282SShri Abhyankar 
1882b4ed282SShri Abhyankar   if (it && !*reason) {
1897fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
1909566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g (relative tolerance)\n", (double)fnorm, (double)snes->ttol));
1912b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
1922b4ed282SShri Abhyankar     }
1932b4ed282SShri Abhyankar   }
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1952b4ed282SShri Abhyankar }
1962b4ed282SShri Abhyankar 
197c1a5e521SShri Abhyankar /*
198c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
199c1a5e521SShri Abhyankar 
200c1a5e521SShri Abhyankar    Input Parameters:
201c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
202c1a5e521SShri Abhyankar 
203c1a5e521SShri Abhyankar    Output Parameters:
204c1a5e521SShri Abhyankar .  X - Bound projected X
205c1a5e521SShri Abhyankar 
206c1a5e521SShri Abhyankar */
207c1a5e521SShri Abhyankar 
208d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X)
209d71ae5a4SJacob Faibussowitsch {
210c1a5e521SShri Abhyankar   const PetscScalar *xl, *xu;
211c1a5e521SShri Abhyankar   PetscScalar       *x;
212c1a5e521SShri Abhyankar   PetscInt           i, n;
213c1a5e521SShri Abhyankar 
214c1a5e521SShri Abhyankar   PetscFunctionBegin;
2159566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
2169566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
2179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
2189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
219c1a5e521SShri Abhyankar 
220c1a5e521SShri Abhyankar   for (i = 0; i < n; i++) {
22110f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
22210f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
223c1a5e521SShri Abhyankar   }
2249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
2259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
2269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
2273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
228c1a5e521SShri Abhyankar }
229c1a5e521SShri Abhyankar 
230ceaaa498SBarry Smith /*@
231e4094ef1SJacob Faibussowitsch   SNESVIGetActiveSetIS - Gets the global indices for the active set variables
232693ddf92SShri Abhyankar 
233ceaaa498SBarry Smith   Input Parameters:
234ceaaa498SBarry Smith + snes - the `SNES` context
235ceaaa498SBarry Smith . X    - the `snes` solution vector
236e4094ef1SJacob Faibussowitsch - F    - the nonlinear function vector
237693ddf92SShri Abhyankar 
238e4094ef1SJacob Faibussowitsch   Output Parameter:
239693ddf92SShri Abhyankar . ISact - active set index set
240ceaaa498SBarry Smith 
241ceaaa498SBarry Smith   Level: developer
242ceaaa498SBarry Smith 
243ceaaa498SBarry Smith .seealso: `SNES`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`
244ceaaa498SBarry Smith @*/
245d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact)
246d71ae5a4SJacob Faibussowitsch {
247c2fc9fa9SBarry Smith   Vec                Xl = snes->xl, Xu = snes->xu;
248693ddf92SShri Abhyankar   const PetscScalar *x, *f, *xl, *xu;
249693ddf92SShri Abhyankar   PetscInt          *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0;
250349187a7SBarry Smith   PetscReal          zerotolerance = snes->vizerotolerance;
251d950fb63SShri Abhyankar 
252d950fb63SShri Abhyankar   PetscFunctionBegin;
2539566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nlocal));
2549566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh));
2559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xl, &xl));
2579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xu, &xu));
2589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
259693ddf92SShri Abhyankar   /* Compute active set size */
260d950fb63SShri Abhyankar   for (i = 0; i < nlocal; i++) {
261349187a7SBarry 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++;
262d950fb63SShri Abhyankar   }
263693ddf92SShri Abhyankar 
2649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc_isact, &idx_act));
265d950fb63SShri Abhyankar 
266693ddf92SShri Abhyankar   /* Set active set indices */
267d950fb63SShri Abhyankar   for (i = 0; i < nlocal; i++) {
268349187a7SBarry 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;
269d950fb63SShri Abhyankar   }
270d950fb63SShri Abhyankar 
271693ddf92SShri Abhyankar   /* Create active set IS */
2729566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact));
273d950fb63SShri Abhyankar 
2749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xl, &xl));
2769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xu, &xu));
2779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
2783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
279d950fb63SShri Abhyankar }
280d950fb63SShri Abhyankar 
281d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm)
282d71ae5a4SJacob Faibussowitsch {
283fe835674SShri Abhyankar   const PetscScalar *x, *xl, *xu, *f;
284fe835674SShri Abhyankar   PetscInt           i, n;
285feb237baSPierre Jolivet   PetscReal          rnorm, zerotolerance = snes->vizerotolerance;
286fe835674SShri Abhyankar 
287fe835674SShri Abhyankar   PetscFunctionBegin;
2889566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
2899566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
2909566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
2919566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2929566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
293fe835674SShri Abhyankar   rnorm = 0.0;
294fe835674SShri Abhyankar   for (i = 0; i < n; 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))) rnorm += PetscRealPart(PetscConj(f[i]) * f[i]);
2968f918527SKarl Rupp   }
2979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
2989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
2999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
3009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
3011c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
30262d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
3033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
304fe835674SShri Abhyankar }
305fe835674SShri Abhyankar 
306ba38deedSJacob Faibussowitsch static PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu)
307d71ae5a4SJacob Faibussowitsch {
30808da532bSDmitry Karpeev   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31108da532bSDmitry Karpeev }
31208da532bSDmitry Karpeev 
3132b4ed282SShri Abhyankar /*
314c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3152b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3162b4ed282SShri Abhyankar 
3172b4ed282SShri Abhyankar    Input Parameter:
3182b4ed282SShri Abhyankar .  snes - the SNES context
3192b4ed282SShri Abhyankar 
3202b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
3212b4ed282SShri Abhyankar 
3222b4ed282SShri Abhyankar    Notes:
3232b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
3242b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
3252b4ed282SShri Abhyankar    the call to SNESSolve().
3262b4ed282SShri Abhyankar  */
327d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_VI(SNES snes)
328d71ae5a4SJacob Faibussowitsch {
3292b4ed282SShri Abhyankar   PetscInt i_start[3], i_end[3];
3302b4ed282SShri Abhyankar 
3312b4ed282SShri Abhyankar   PetscFunctionBegin;
3329566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 1));
3339566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
3342b4ed282SShri Abhyankar 
33508da532bSDmitry Karpeev   if (!snes->ops->computevariablebounds && snes->dm) {
336a201590fSDmitry Karpeev     PetscBool flag;
3379566063dSJacob Faibussowitsch     PetscCall(DMHasVariableBounds(snes->dm, &flag));
338ad540459SPierre Jolivet     if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
33936b2a9f3SBarry Smith   }
340a201590fSDmitry Karpeev   if (!snes->usersetbounds) {
341c2fc9fa9SBarry Smith     if (snes->ops->computevariablebounds) {
3429566063dSJacob Faibussowitsch       if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3439566063dSJacob Faibussowitsch       if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
344dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
3451aa26658SKarl Rupp     } else if (!snes->xl && !snes->xu) {
3462176524cSBarry Smith       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
3479566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3489566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xl, PETSC_NINFINITY));
3499566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
3509566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xu, PETSC_INFINITY));
3512b4ed282SShri Abhyankar     } else {
3522b4ed282SShri Abhyankar       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
3539566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end));
3549566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1));
3559566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2));
3562b4ed282SShri 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]))
3572b4ed282SShri 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.");
3582b4ed282SShri Abhyankar     }
359a201590fSDmitry Karpeev   }
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3612b4ed282SShri Abhyankar }
362d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_VI(SNES snes)
363d71ae5a4SJacob Faibussowitsch {
3642176524cSBarry Smith   PetscFunctionBegin;
3659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
3669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
3672d6615e8SDmitry Karpeev   snes->usersetbounds = PETSC_FALSE;
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3692176524cSBarry Smith }
3702176524cSBarry Smith 
3712b4ed282SShri Abhyankar /*
3722b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
3732b4ed282SShri Abhyankar    with SNESCreate_VI().
3742b4ed282SShri Abhyankar 
3752b4ed282SShri Abhyankar    Input Parameter:
3762b4ed282SShri Abhyankar .  snes - the SNES context
3772b4ed282SShri Abhyankar 
3782b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
3792b4ed282SShri Abhyankar  */
380d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_VI(SNES snes)
381d71ae5a4SJacob Faibussowitsch {
3822b4ed282SShri Abhyankar   PetscFunctionBegin;
3839566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3842b4ed282SShri Abhyankar 
3852b4ed282SShri Abhyankar   /* clear composed functions */
3862e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL));
3872e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL));
3883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3892b4ed282SShri Abhyankar }
3907fe79bc4SShri Abhyankar 
3915559b345SBarry Smith /*@
392f0b84518SBarry Smith   SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving
393f6dfbefdSBarry Smith   (differential) variable inequalities.
3942b4ed282SShri Abhyankar 
3952b4ed282SShri Abhyankar   Input Parameters:
396f6dfbefdSBarry Smith + snes - the `SNES` context.
3972b4ed282SShri Abhyankar . xl   - lower bound.
398a2b725a8SWilliam Gropp - xu   - upper bound.
3992b4ed282SShri Abhyankar 
4002fe279fdSBarry Smith   Level: advanced
4012fe279fdSBarry Smith 
4022b4ed282SShri Abhyankar   Notes:
4032b4ed282SShri Abhyankar   If this routine is not called then the lower and upper bounds are set to
404f6dfbefdSBarry Smith   `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`.
40584c105d7SBarry Smith 
406f0b84518SBarry Smith   Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.
407f0b84518SBarry Smith 
408f6dfbefdSBarry Smith   For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY`
409f6dfbefdSBarry Smith 
41033e0caf2SBarry Smith   `SNESVISetComputeVariableBounds()` can be used to provide a function that computes the bounds. This should be used if you are using, for example, grid
411f6dfbefdSBarry Smith   sequencing and need bounds set for a variety of vectors
412f6dfbefdSBarry Smith 
413*945b2ebdSBarry Smith .seealso: [](sec_vi), `SNES`, `SNESVIGetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `'SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
4145559b345SBarry Smith @*/
415d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
416d71ae5a4SJacob Faibussowitsch {
4175f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES, Vec, Vec);
4182b4ed282SShri Abhyankar 
4192b4ed282SShri Abhyankar   PetscFunctionBegin;
4202b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
4212b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl, VEC_CLASSID, 2);
4222b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu, VEC_CLASSID, 3);
4239566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f));
424cac4c232SBarry Smith   if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu));
4259566063dSJacob Faibussowitsch   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
426a201590fSDmitry Karpeev   snes->usersetbounds = PETSC_TRUE;
4273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42861589011SJed Brown }
42961589011SJed Brown 
430d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu)
431d71ae5a4SJacob Faibussowitsch {
43261589011SJed Brown   const PetscScalar *xxl, *xxu;
43361589011SJed Brown   PetscInt           i, n, cnt = 0;
43461589011SJed Brown 
43561589011SJed Brown   PetscFunctionBegin;
4369566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL));
4375f80ce2aSJacob Faibussowitsch   PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first");
438077aedafSJed Brown   {
439077aedafSJed Brown     PetscInt xlN, xuN, N;
4409566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xl, &xlN));
4419566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xu, &xuN));
4429566063dSJacob Faibussowitsch     PetscCall(VecGetSize(snes->vec_func, &N));
44363a3b9bcSJacob Faibussowitsch     PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N);
44463a3b9bcSJacob Faibussowitsch     PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N);
445077aedafSJed Brown   }
4469566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xl));
4479566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xu));
4489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
4499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
450c2fc9fa9SBarry Smith   snes->xl = xl;
451c2fc9fa9SBarry Smith   snes->xu = xu;
4529566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xl, &n));
4539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xl, &xxl));
4549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xu, &xxu));
455e270355aSBarry Smith   for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
4561aa26658SKarl Rupp 
4571c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
4589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xl, &xxl));
4599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xu, &xxu));
4603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4612b4ed282SShri Abhyankar }
46292c02d66SPeter Brune 
463cf836535SBlaise Bourdin /*@
464*945b2ebdSBarry Smith   SNESVIGetVariableBounds - Gets the lower and upper bounds for the solution vector. xl <= x <= xu. These are used in solving
465cf836535SBlaise Bourdin   (differential) variable inequalities.
466cf836535SBlaise Bourdin 
467cf836535SBlaise Bourdin   Input Parameters:
468cf836535SBlaise Bourdin + snes - the `SNES` context.
469cf836535SBlaise Bourdin . xl   - lower bound (may be `NULL`)
470cf836535SBlaise Bourdin - xu   - upper bound (may be `NULL`)
471cf836535SBlaise Bourdin 
4722fe279fdSBarry Smith   Level: advanced
4732fe279fdSBarry Smith 
474cf836535SBlaise Bourdin   Notes:
475cf836535SBlaise Bourdin   These vectors are owned by the `SNESVI` and should not be destroyed by the caller
476cf836535SBlaise Bourdin 
477*945b2ebdSBarry Smith .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `'SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
478cf836535SBlaise Bourdin @*/
479cf836535SBlaise Bourdin PetscErrorCode SNESVIGetVariableBounds(SNES snes, Vec *xl, Vec *xu)
480cf836535SBlaise Bourdin {
481cf836535SBlaise Bourdin   PetscFunctionBegin;
482cf836535SBlaise Bourdin   PetscCheck(snes->usersetbounds, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must set SNESVI bounds before calling SNESVIGetVariableBounds()");
483cf836535SBlaise Bourdin   if (xl) *xl = snes->xl;
484cf836535SBlaise Bourdin   if (xu) *xu = snes->xu;
4853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
486cf836535SBlaise Bourdin }
487cf836535SBlaise Bourdin 
488d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject)
489d71ae5a4SJacob Faibussowitsch {
4908afaa268SBarry Smith   PetscBool flg = PETSC_FALSE;
4912b4ed282SShri Abhyankar 
4922b4ed282SShri Abhyankar   PetscFunctionBegin;
493d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options");
4949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL));
4959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL));
4961baa6e33SBarry Smith   if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL));
497de34d3e9SBarry Smith   flg = PETSC_FALSE;
4989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor_residual", "Monitor residual all non-active variables; using zero for active constraints", "SNESMonitorVIResidual", flg, &flg, NULL));
4991baa6e33SBarry Smith   if (flg) PetscCall(SNESMonitorSet(snes, SNESVIMonitorResidual, PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)), NULL));
500d0609cedSBarry Smith   PetscOptionsHeadEnd();
5013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
502ed70e96aSJungho Lee }
503