xref: /petsc/src/snes/impls/vi/vi.c (revision b6f515db33ea3628673fe08d569f5d88d953eda2)
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 @*/
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 
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 
51*b6f515dbSMatthew 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;
57*b6f515dbSMatthew 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));
62*b6f515dbSMatthew 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));
66*b6f515dbSMatthew G. Knepley   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
67*b6f515dbSMatthew G. Knepley   PetscCall(VecView(Finactive, vf->viewer));
68*b6f515dbSMatthew G. Knepley   PetscCall(PetscViewerPopFormat(vf->viewer));
69*b6f515dbSMatthew G. Knepley   PetscCall(PetscObjectCompose((PetscObject)Finactive, "__Vec_bc_zero__", NULL));
709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Finactive));
713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72ffdf2a20SBarry Smith }
73ffdf2a20SBarry Smith 
74*b6f515dbSMatthew G. Knepley static PetscErrorCode SNESVIMonitorActive(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
75*b6f515dbSMatthew G. Knepley {
76*b6f515dbSMatthew G. Knepley   Vec X, F, A;
77*b6f515dbSMatthew G. Knepley   IS  isactive;
78*b6f515dbSMatthew G. Knepley 
79*b6f515dbSMatthew G. Knepley   PetscFunctionBegin;
80*b6f515dbSMatthew G. Knepley   PetscValidHeaderSpecific(vf->viewer, PETSC_VIEWER_CLASSID, 4);
81*b6f515dbSMatthew G. Knepley   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
82*b6f515dbSMatthew G. Knepley   PetscCall(SNESGetSolution(snes, &X));
83*b6f515dbSMatthew G. Knepley   PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive));
84*b6f515dbSMatthew G. Knepley   PetscCall(VecDuplicate(F, &A));
85*b6f515dbSMatthew G. Knepley   PetscCall(PetscObjectCompose((PetscObject)A, "__Vec_bc_zero__", (PetscObject)snes));
86*b6f515dbSMatthew G. Knepley   PetscCall(VecSet(A, 0.));
87*b6f515dbSMatthew G. Knepley   PetscCall(VecISSet(A, isactive, 1.));
88*b6f515dbSMatthew G. Knepley   PetscCall(ISDestroy(&isactive));
89*b6f515dbSMatthew G. Knepley   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
90*b6f515dbSMatthew G. Knepley   PetscCall(VecView(A, vf->viewer));
91*b6f515dbSMatthew G. Knepley   PetscCall(PetscViewerPopFormat(vf->viewer));
92*b6f515dbSMatthew G. Knepley   PetscCall(PetscObjectCompose((PetscObject)A, "__Vec_bc_zero__", NULL));
93*b6f515dbSMatthew G. Knepley   PetscCall(VecDestroy(&A));
94*b6f515dbSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
95*b6f515dbSMatthew G. Knepley }
96*b6f515dbSMatthew G. Knepley 
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 */
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 */
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"));
2062b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
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 
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 @*/
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 
308d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm)
309d71ae5a4SJacob Faibussowitsch {
310fe835674SShri Abhyankar   const PetscScalar *x, *xl, *xu, *f;
311fe835674SShri Abhyankar   PetscInt           i, n;
312feb237baSPierre Jolivet   PetscReal          rnorm, zerotolerance = snes->vizerotolerance;
313fe835674SShri Abhyankar 
314fe835674SShri Abhyankar   PetscFunctionBegin;
3159566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
3169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
3179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
3189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
3199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
320fe835674SShri Abhyankar   rnorm = 0.0;
321fe835674SShri Abhyankar   for (i = 0; i < n; i++) {
322f4f49eeaSPierre 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]);
3238f918527SKarl Rupp   }
3249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
3259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
3269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
3279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
328462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
32962d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
3303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
331fe835674SShri Abhyankar }
332fe835674SShri Abhyankar 
333ba38deedSJacob Faibussowitsch static PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu)
334d71ae5a4SJacob Faibussowitsch {
33508da532bSDmitry Karpeev   PetscFunctionBegin;
3369566063dSJacob Faibussowitsch   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33808da532bSDmitry Karpeev }
33908da532bSDmitry Karpeev 
3402b4ed282SShri Abhyankar /*
341c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3422b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3432b4ed282SShri Abhyankar 
3442b4ed282SShri Abhyankar    Input Parameter:
3452b4ed282SShri Abhyankar .  snes - the SNES context
3462b4ed282SShri Abhyankar 
3472b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
3482b4ed282SShri Abhyankar 
3492b4ed282SShri Abhyankar    Notes:
3502b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
3512b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
3522b4ed282SShri Abhyankar    the call to SNESSolve().
3532b4ed282SShri Abhyankar  */
354d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_VI(SNES snes)
355d71ae5a4SJacob Faibussowitsch {
3562b4ed282SShri Abhyankar   PetscInt i_start[3], i_end[3];
3572b4ed282SShri Abhyankar 
3582b4ed282SShri Abhyankar   PetscFunctionBegin;
3599566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 1));
3609566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
3612b4ed282SShri Abhyankar 
36208da532bSDmitry Karpeev   if (!snes->ops->computevariablebounds && snes->dm) {
363a201590fSDmitry Karpeev     PetscBool flag;
3649566063dSJacob Faibussowitsch     PetscCall(DMHasVariableBounds(snes->dm, &flag));
365ad540459SPierre Jolivet     if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
36636b2a9f3SBarry Smith   }
367a201590fSDmitry Karpeev   if (!snes->usersetbounds) {
368c2fc9fa9SBarry Smith     if (snes->ops->computevariablebounds) {
3699566063dSJacob Faibussowitsch       if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3709566063dSJacob Faibussowitsch       if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
371dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
3721aa26658SKarl Rupp     } else if (!snes->xl && !snes->xu) {
3732176524cSBarry Smith       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
3749566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3759566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xl, PETSC_NINFINITY));
3769566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
3779566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xu, PETSC_INFINITY));
3782b4ed282SShri Abhyankar     } else {
3792b4ed282SShri Abhyankar       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
3809566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end));
3819566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1));
3829566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2));
3832b4ed282SShri 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]))
3842b4ed282SShri 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.");
3852b4ed282SShri Abhyankar     }
386a201590fSDmitry Karpeev   }
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3882b4ed282SShri Abhyankar }
389d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_VI(SNES snes)
390d71ae5a4SJacob Faibussowitsch {
3912176524cSBarry Smith   PetscFunctionBegin;
3929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
3939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
3942d6615e8SDmitry Karpeev   snes->usersetbounds = PETSC_FALSE;
3953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3962176524cSBarry Smith }
3972176524cSBarry Smith 
3982b4ed282SShri Abhyankar /*
3992b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
4002b4ed282SShri Abhyankar    with SNESCreate_VI().
4012b4ed282SShri Abhyankar 
4022b4ed282SShri Abhyankar    Input Parameter:
4032b4ed282SShri Abhyankar .  snes - the SNES context
4042b4ed282SShri Abhyankar 
4052b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4062b4ed282SShri Abhyankar  */
407d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_VI(SNES snes)
408d71ae5a4SJacob Faibussowitsch {
4092b4ed282SShri Abhyankar   PetscFunctionBegin;
4109566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
4112b4ed282SShri Abhyankar 
4122b4ed282SShri Abhyankar   /* clear composed functions */
4132e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL));
4142e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL));
4153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4162b4ed282SShri Abhyankar }
4177fe79bc4SShri Abhyankar 
4185559b345SBarry Smith /*@
419420bcc1bSBarry Smith   SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. `xl` <= x <= `xu`. This allows solving
420f6dfbefdSBarry Smith   (differential) variable inequalities.
4212b4ed282SShri Abhyankar 
4222b4ed282SShri Abhyankar   Input Parameters:
423f6dfbefdSBarry Smith + snes - the `SNES` context.
4242b4ed282SShri Abhyankar . xl   - lower bound.
425a2b725a8SWilliam Gropp - xu   - upper bound.
4262b4ed282SShri Abhyankar 
4272fe279fdSBarry Smith   Level: advanced
4282fe279fdSBarry Smith 
4292b4ed282SShri Abhyankar   Notes:
4302b4ed282SShri Abhyankar   If this routine is not called then the lower and upper bounds are set to
431f6dfbefdSBarry Smith   `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`.
43284c105d7SBarry Smith 
433420bcc1bSBarry Smith   Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS` or semi-smooth `SNESVINEWTONSSLS` solvers.
434f0b84518SBarry Smith 
435f6dfbefdSBarry Smith   For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY`
436f6dfbefdSBarry Smith 
43733e0caf2SBarry Smith   `SNESVISetComputeVariableBounds()` can be used to provide a function that computes the bounds. This should be used if you are using, for example, grid
438f6dfbefdSBarry Smith   sequencing and need bounds set for a variety of vectors
439f6dfbefdSBarry Smith 
4402913281cSPierre Jolivet .seealso: [](sec_vi), `SNES`, `SNESVIGetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
4415559b345SBarry Smith @*/
442d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
443d71ae5a4SJacob Faibussowitsch {
4445f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES, Vec, Vec);
4452b4ed282SShri Abhyankar 
4462b4ed282SShri Abhyankar   PetscFunctionBegin;
4472b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
4482b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl, VEC_CLASSID, 2);
4492b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu, VEC_CLASSID, 3);
4509566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f));
451cac4c232SBarry Smith   if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu));
4529566063dSJacob Faibussowitsch   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
453a201590fSDmitry Karpeev   snes->usersetbounds = PETSC_TRUE;
4543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45561589011SJed Brown }
45661589011SJed Brown 
457d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu)
458d71ae5a4SJacob Faibussowitsch {
45961589011SJed Brown   const PetscScalar *xxl, *xxu;
46061589011SJed Brown   PetscInt           i, n, cnt = 0;
46161589011SJed Brown 
46261589011SJed Brown   PetscFunctionBegin;
4639566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL));
4645f80ce2aSJacob Faibussowitsch   PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first");
465077aedafSJed Brown   {
466077aedafSJed Brown     PetscInt xlN, xuN, N;
4679566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xl, &xlN));
4689566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xu, &xuN));
4699566063dSJacob Faibussowitsch     PetscCall(VecGetSize(snes->vec_func, &N));
47063a3b9bcSJacob Faibussowitsch     PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N);
47163a3b9bcSJacob Faibussowitsch     PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N);
472077aedafSJed Brown   }
4739566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xl));
4749566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xu));
4759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
4769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
477c2fc9fa9SBarry Smith   snes->xl = xl;
478c2fc9fa9SBarry Smith   snes->xu = xu;
4799566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xl, &n));
4809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xl, &xxl));
4819566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xu, &xxu));
482e270355aSBarry Smith   for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
4831aa26658SKarl Rupp 
484462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
4859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xl, &xxl));
4869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xu, &xxu));
4873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4882b4ed282SShri Abhyankar }
48992c02d66SPeter Brune 
490cf836535SBlaise Bourdin /*@
491420bcc1bSBarry Smith   SNESVIGetVariableBounds - Gets the lower and upper bounds for the solution vector. `xl` <= x <= `xu`. These are used in solving
492cf836535SBlaise Bourdin   (differential) variable inequalities.
493cf836535SBlaise Bourdin 
494cf836535SBlaise Bourdin   Input Parameters:
495cf836535SBlaise Bourdin + snes - the `SNES` context.
496cf836535SBlaise Bourdin . xl   - lower bound (may be `NULL`)
497cf836535SBlaise Bourdin - xu   - upper bound (may be `NULL`)
498cf836535SBlaise Bourdin 
4992fe279fdSBarry Smith   Level: advanced
5002fe279fdSBarry Smith 
501420bcc1bSBarry Smith   Note:
502cf836535SBlaise Bourdin   These vectors are owned by the `SNESVI` and should not be destroyed by the caller
503cf836535SBlaise Bourdin 
5042913281cSPierre Jolivet .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY`
505cf836535SBlaise Bourdin @*/
506cf836535SBlaise Bourdin PetscErrorCode SNESVIGetVariableBounds(SNES snes, Vec *xl, Vec *xu)
507cf836535SBlaise Bourdin {
508cf836535SBlaise Bourdin   PetscFunctionBegin;
509cf836535SBlaise Bourdin   PetscCheck(snes->usersetbounds, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must set SNESVI bounds before calling SNESVIGetVariableBounds()");
510cf836535SBlaise Bourdin   if (xl) *xl = snes->xl;
511cf836535SBlaise Bourdin   if (xu) *xu = snes->xu;
5123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
513cf836535SBlaise Bourdin }
514cf836535SBlaise Bourdin 
515d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject)
516d71ae5a4SJacob Faibussowitsch {
5178afaa268SBarry Smith   PetscBool flg = PETSC_FALSE;
5182b4ed282SShri Abhyankar 
5192b4ed282SShri Abhyankar   PetscFunctionBegin;
520d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options");
5219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL));
5229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL));
5231baa6e33SBarry Smith   if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL));
524de34d3e9SBarry Smith   flg = PETSC_FALSE;
525*b6f515dbSMatthew G. Knepley   PetscCall(SNESMonitorSetFromOptions(snes, "-snes_vi_monitor_residual", "View residual at each iteration, using zero for active constraints", "SNESVIMonitorResidual", SNESVIMonitorResidual, NULL));
526*b6f515dbSMatthew G. Knepley   PetscCall(SNESMonitorSetFromOptions(snes, "-snes_vi_monitor_active", "View active set at each iteration, using zero for inactive dofs", "SNESVIMonitorActive", SNESVIMonitorActive, NULL));
527d0609cedSBarry Smith   PetscOptionsHeadEnd();
5283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
529ed70e96aSJungho Lee }
530