xref: /petsc/src/snes/impls/vi/vi.c (revision ceaaa4989964adb3f5eb130cb04b8f49c83e49c9)
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`,
30e4094ef1SJacob Faibussowitsch           `'SNESSetType()`
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 
51d71ae5a4SJacob Faibussowitsch 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 
70d71ae5a4SJacob Faibussowitsch 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 /*
1582b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
1592b4ed282SShri Abhyankar */
160d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVICheckResidual_Private(SNES snes, Mat A, Vec F, Vec X, Vec W1, Vec W2)
161d71ae5a4SJacob Faibussowitsch {
1622b4ed282SShri Abhyankar   PetscReal a1, a2;
163ace3abfcSBarry Smith   PetscBool hastranspose;
1642b4ed282SShri Abhyankar 
1652b4ed282SShri Abhyankar   PetscFunctionBegin;
1669566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
1672b4ed282SShri Abhyankar   if (hastranspose) {
1689566063dSJacob Faibussowitsch     PetscCall(MatMult(A, X, W1));
1699566063dSJacob Faibussowitsch     PetscCall(VecAXPY(W1, -1.0, F));
1702b4ed282SShri Abhyankar 
1712b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1729566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, W1, W2));
1739566063dSJacob Faibussowitsch     PetscCall(VecNorm(W1, NORM_2, &a1));
1749566063dSJacob Faibussowitsch     PetscCall(VecNorm(W2, NORM_2, &a2));
17548a46eb9SPierre Jolivet     if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n", (double)(a2 / a1)));
1762b4ed282SShri Abhyankar   }
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1782b4ed282SShri Abhyankar }
1792b4ed282SShri Abhyankar 
1802b4ed282SShri Abhyankar /*
1818d359177SBarry Smith   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
1822b4ed282SShri Abhyankar 
1832b4ed282SShri Abhyankar   Notes:
1842b4ed282SShri Abhyankar   The convergence criterion currently implemented is
1852b4ed282SShri Abhyankar   merit < abstol
1862b4ed282SShri Abhyankar   merit < rtol*merit_initial
1872b4ed282SShri Abhyankar */
188d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESConvergedDefault_VI(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gradnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
189d71ae5a4SJacob Faibussowitsch {
1902b4ed282SShri Abhyankar   PetscFunctionBegin;
1912b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1924f572ea9SToby Isaac   PetscAssertPointer(reason, 6);
1932b4ed282SShri Abhyankar 
1942b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
1952b4ed282SShri Abhyankar 
1962b4ed282SShri Abhyankar   if (!it) {
1972b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
1987fe79bc4SShri Abhyankar     snes->ttol = fnorm * snes->rtol;
1992b4ed282SShri Abhyankar   }
2007fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
2019566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n"));
2022b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
2030d6f27a8SBarry Smith   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
2049566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g\n", (double)fnorm, (double)snes->abstol));
2052b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
206e71169deSBarry Smith   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
20763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs));
2082b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
2092b4ed282SShri Abhyankar   }
2102b4ed282SShri Abhyankar 
2112b4ed282SShri Abhyankar   if (it && !*reason) {
2127fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2139566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g (relative tolerance)\n", (double)fnorm, (double)snes->ttol));
2142b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2152b4ed282SShri Abhyankar     }
2162b4ed282SShri Abhyankar   }
2173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2182b4ed282SShri Abhyankar }
2192b4ed282SShri Abhyankar 
220c1a5e521SShri Abhyankar /*
221c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
222c1a5e521SShri Abhyankar 
223c1a5e521SShri Abhyankar    Input Parameters:
224c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
225c1a5e521SShri Abhyankar 
226c1a5e521SShri Abhyankar    Output Parameters:
227c1a5e521SShri Abhyankar .  X - Bound projected X
228c1a5e521SShri Abhyankar 
229c1a5e521SShri Abhyankar */
230c1a5e521SShri Abhyankar 
231d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X)
232d71ae5a4SJacob Faibussowitsch {
233c1a5e521SShri Abhyankar   const PetscScalar *xl, *xu;
234c1a5e521SShri Abhyankar   PetscScalar       *x;
235c1a5e521SShri Abhyankar   PetscInt           i, n;
236c1a5e521SShri Abhyankar 
237c1a5e521SShri Abhyankar   PetscFunctionBegin;
2389566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
2399566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
2409566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
2419566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
242c1a5e521SShri Abhyankar 
243c1a5e521SShri Abhyankar   for (i = 0; i < n; i++) {
24410f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
24510f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
246c1a5e521SShri Abhyankar   }
2479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
2489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
2499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
2503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
251c1a5e521SShri Abhyankar }
252c1a5e521SShri Abhyankar 
253*ceaaa498SBarry Smith /*@
254e4094ef1SJacob Faibussowitsch   SNESVIGetActiveSetIS - Gets the global indices for the active set variables
255693ddf92SShri Abhyankar 
256*ceaaa498SBarry Smith   Input Parameters:
257*ceaaa498SBarry Smith + snes - the `SNES` context
258*ceaaa498SBarry Smith . X    - the `snes` solution vector
259e4094ef1SJacob Faibussowitsch - F    - the nonlinear function vector
260693ddf92SShri Abhyankar 
261e4094ef1SJacob Faibussowitsch   Output Parameter:
262693ddf92SShri Abhyankar . ISact - active set index set
263*ceaaa498SBarry Smith 
264*ceaaa498SBarry Smith   Level: developer
265*ceaaa498SBarry Smith 
266*ceaaa498SBarry Smith .seealso: `SNES`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`
267*ceaaa498SBarry Smith @*/
268d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact)
269d71ae5a4SJacob Faibussowitsch {
270c2fc9fa9SBarry Smith   Vec                Xl = snes->xl, Xu = snes->xu;
271693ddf92SShri Abhyankar   const PetscScalar *x, *f, *xl, *xu;
272693ddf92SShri Abhyankar   PetscInt          *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0;
273349187a7SBarry Smith   PetscReal          zerotolerance = snes->vizerotolerance;
274d950fb63SShri Abhyankar 
275d950fb63SShri Abhyankar   PetscFunctionBegin;
2769566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &nlocal));
2779566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh));
2789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2799566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xl, &xl));
2809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xu, &xu));
2819566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
282693ddf92SShri Abhyankar   /* Compute active set size */
283d950fb63SShri Abhyankar   for (i = 0; i < nlocal; i++) {
284349187a7SBarry 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++;
285d950fb63SShri Abhyankar   }
286693ddf92SShri Abhyankar 
2879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc_isact, &idx_act));
288d950fb63SShri Abhyankar 
289693ddf92SShri Abhyankar   /* Set active set indices */
290d950fb63SShri Abhyankar   for (i = 0; i < nlocal; i++) {
291349187a7SBarry 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;
292d950fb63SShri Abhyankar   }
293d950fb63SShri Abhyankar 
294693ddf92SShri Abhyankar   /* Create active set IS */
2959566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact));
296d950fb63SShri Abhyankar 
2979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xl, &xl));
2999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xu, &xu));
3009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
3013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
302d950fb63SShri Abhyankar }
303d950fb63SShri Abhyankar 
304d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVICreateIndexSets_RS(SNES snes, Vec X, Vec F, IS *ISact, IS *ISinact)
305d71ae5a4SJacob Faibussowitsch {
306077aedafSJed Brown   PetscInt rstart, rend;
307693ddf92SShri Abhyankar 
308693ddf92SShri Abhyankar   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(SNESVIGetActiveSetIS(snes, X, F, ISact));
3109566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &rstart, &rend));
3119566063dSJacob Faibussowitsch   PetscCall(ISComplement(*ISact, rstart, rend, ISinact));
3123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
313693ddf92SShri Abhyankar }
314693ddf92SShri Abhyankar 
315d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm)
316d71ae5a4SJacob Faibussowitsch {
317fe835674SShri Abhyankar   const PetscScalar *x, *xl, *xu, *f;
318fe835674SShri Abhyankar   PetscInt           i, n;
319feb237baSPierre Jolivet   PetscReal          rnorm, zerotolerance = snes->vizerotolerance;
320fe835674SShri Abhyankar 
321fe835674SShri Abhyankar   PetscFunctionBegin;
3229566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
3239566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xl, &xl));
3249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(snes->xu, &xu));
3259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
3269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &f));
327fe835674SShri Abhyankar   rnorm = 0.0;
328fe835674SShri Abhyankar   for (i = 0; i < n; i++) {
329349187a7SBarry 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]);
3308f918527SKarl Rupp   }
3319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &f));
3329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xl, &xl));
3339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(snes->xu, &xu));
3349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
3351c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
33662d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
338fe835674SShri Abhyankar }
339fe835674SShri Abhyankar 
340d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu)
341d71ae5a4SJacob Faibussowitsch {
34208da532bSDmitry Karpeev   PetscFunctionBegin;
3439566063dSJacob Faibussowitsch   PetscCall(DMComputeVariableBounds(snes->dm, xl, xu));
3443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34508da532bSDmitry Karpeev }
34608da532bSDmitry Karpeev 
3472b4ed282SShri Abhyankar /*
348c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3492b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3502b4ed282SShri Abhyankar 
3512b4ed282SShri Abhyankar    Input Parameter:
3522b4ed282SShri Abhyankar .  snes - the SNES context
3532b4ed282SShri Abhyankar 
3542b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
3552b4ed282SShri Abhyankar 
3562b4ed282SShri Abhyankar    Notes:
3572b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
3582b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
3592b4ed282SShri Abhyankar    the call to SNESSolve().
3602b4ed282SShri Abhyankar  */
361d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_VI(SNES snes)
362d71ae5a4SJacob Faibussowitsch {
3632b4ed282SShri Abhyankar   PetscInt i_start[3], i_end[3];
3642b4ed282SShri Abhyankar 
3652b4ed282SShri Abhyankar   PetscFunctionBegin;
3669566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 1));
3679566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
3682b4ed282SShri Abhyankar 
36908da532bSDmitry Karpeev   if (!snes->ops->computevariablebounds && snes->dm) {
370a201590fSDmitry Karpeev     PetscBool flag;
3719566063dSJacob Faibussowitsch     PetscCall(DMHasVariableBounds(snes->dm, &flag));
372ad540459SPierre Jolivet     if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
37336b2a9f3SBarry Smith   }
374a201590fSDmitry Karpeev   if (!snes->usersetbounds) {
375c2fc9fa9SBarry Smith     if (snes->ops->computevariablebounds) {
3769566063dSJacob Faibussowitsch       if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3779566063dSJacob Faibussowitsch       if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
378dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu);
3791aa26658SKarl Rupp     } else if (!snes->xl && !snes->xu) {
3802176524cSBarry Smith       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
3819566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xl));
3829566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xl, PETSC_NINFINITY));
3839566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol, &snes->xu));
3849566063dSJacob Faibussowitsch       PetscCall(VecSet(snes->xu, PETSC_INFINITY));
3852b4ed282SShri Abhyankar     } else {
3862b4ed282SShri Abhyankar       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
3879566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end));
3889566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1));
3899566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2));
3902b4ed282SShri 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]))
3912b4ed282SShri 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.");
3922b4ed282SShri Abhyankar     }
393a201590fSDmitry Karpeev   }
3943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3952b4ed282SShri Abhyankar }
396d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_VI(SNES snes)
397d71ae5a4SJacob Faibussowitsch {
3982176524cSBarry Smith   PetscFunctionBegin;
3999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
4009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
4012d6615e8SDmitry Karpeev   snes->usersetbounds = PETSC_FALSE;
4023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4032176524cSBarry Smith }
4042176524cSBarry Smith 
4052b4ed282SShri Abhyankar /*
4062b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
4072b4ed282SShri Abhyankar    with SNESCreate_VI().
4082b4ed282SShri Abhyankar 
4092b4ed282SShri Abhyankar    Input Parameter:
4102b4ed282SShri Abhyankar .  snes - the SNES context
4112b4ed282SShri Abhyankar 
4122b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4132b4ed282SShri Abhyankar  */
414d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_VI(SNES snes)
415d71ae5a4SJacob Faibussowitsch {
4162b4ed282SShri Abhyankar   PetscFunctionBegin;
4179566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
4182b4ed282SShri Abhyankar 
4192b4ed282SShri Abhyankar   /* clear composed functions */
4202e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL));
4212e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL));
4223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4232b4ed282SShri Abhyankar }
4247fe79bc4SShri Abhyankar 
4255559b345SBarry Smith /*@
426f0b84518SBarry Smith   SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving
427f6dfbefdSBarry Smith   (differential) variable inequalities.
4282b4ed282SShri Abhyankar 
4292b4ed282SShri Abhyankar   Input Parameters:
430f6dfbefdSBarry Smith + snes - the `SNES` context.
4312b4ed282SShri Abhyankar . xl   - lower bound.
432a2b725a8SWilliam Gropp - xu   - upper bound.
4332b4ed282SShri Abhyankar 
4342fe279fdSBarry Smith   Level: advanced
4352fe279fdSBarry Smith 
4362b4ed282SShri Abhyankar   Notes:
4372b4ed282SShri Abhyankar   If this routine is not called then the lower and upper bounds are set to
438f6dfbefdSBarry Smith   `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`.
43984c105d7SBarry Smith 
440f0b84518SBarry Smith   Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers.
441f0b84518SBarry Smith 
442f6dfbefdSBarry Smith   For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY`
443f6dfbefdSBarry Smith 
44433e0caf2SBarry Smith   `SNESVISetComputeVariableBounds()` can be used to provide a function that computes the bounds. This should be used if you are using, for example, grid
445f6dfbefdSBarry Smith   sequencing and need bounds set for a variety of vectors
446f6dfbefdSBarry Smith 
447e4094ef1SJacob Faibussowitsch .seealso: [](sec_vi), `SNES`, `SNESVIGetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `'SNESSetType()`
4485559b345SBarry Smith @*/
449d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
450d71ae5a4SJacob Faibussowitsch {
4515f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(SNES, Vec, Vec);
4522b4ed282SShri Abhyankar 
4532b4ed282SShri Abhyankar   PetscFunctionBegin;
4542b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
4552b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl, VEC_CLASSID, 2);
4562b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu, VEC_CLASSID, 3);
4579566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f));
458cac4c232SBarry Smith   if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu));
4599566063dSJacob Faibussowitsch   else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu));
460a201590fSDmitry Karpeev   snes->usersetbounds = PETSC_TRUE;
4613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46261589011SJed Brown }
46361589011SJed Brown 
464d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu)
465d71ae5a4SJacob Faibussowitsch {
46661589011SJed Brown   const PetscScalar *xxl, *xxu;
46761589011SJed Brown   PetscInt           i, n, cnt = 0;
46861589011SJed Brown 
46961589011SJed Brown   PetscFunctionBegin;
4709566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL));
4715f80ce2aSJacob Faibussowitsch   PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first");
472077aedafSJed Brown   {
473077aedafSJed Brown     PetscInt xlN, xuN, N;
4749566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xl, &xlN));
4759566063dSJacob Faibussowitsch     PetscCall(VecGetSize(xu, &xuN));
4769566063dSJacob Faibussowitsch     PetscCall(VecGetSize(snes->vec_func, &N));
47763a3b9bcSJacob Faibussowitsch     PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N);
47863a3b9bcSJacob Faibussowitsch     PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N);
479077aedafSJed Brown   }
4809566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xl));
4819566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)xu));
4829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xl));
4839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&snes->xu));
484c2fc9fa9SBarry Smith   snes->xl = xl;
485c2fc9fa9SBarry Smith   snes->xu = xu;
4869566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xl, &n));
4879566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xl, &xxl));
4889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xu, &xxu));
489e270355aSBarry Smith   for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
4901aa26658SKarl Rupp 
4911c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes)));
4929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xl, &xxl));
4939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xu, &xxu));
4943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4952b4ed282SShri Abhyankar }
49692c02d66SPeter Brune 
497cf836535SBlaise Bourdin /*@
498cf836535SBlaise Bourdin   SNESVIGetVariableBounds - Gets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving
499cf836535SBlaise Bourdin   (differential) variable inequalities.
500cf836535SBlaise Bourdin 
501cf836535SBlaise Bourdin   Input Parameters:
502cf836535SBlaise Bourdin + snes - the `SNES` context.
503cf836535SBlaise Bourdin . xl   - lower bound (may be `NULL`)
504cf836535SBlaise Bourdin - xu   - upper bound (may be `NULL`)
505cf836535SBlaise Bourdin 
5062fe279fdSBarry Smith   Level: advanced
5072fe279fdSBarry Smith 
508cf836535SBlaise Bourdin   Notes:
509cf836535SBlaise Bourdin   These vectors are owned by the `SNESVI` and should not be destroyed by the caller
510cf836535SBlaise Bourdin 
511e4094ef1SJacob Faibussowitsch .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `'SNESSetType()`
512cf836535SBlaise Bourdin @*/
513cf836535SBlaise Bourdin PetscErrorCode SNESVIGetVariableBounds(SNES snes, Vec *xl, Vec *xu)
514cf836535SBlaise Bourdin {
515cf836535SBlaise Bourdin   PetscFunctionBegin;
516cf836535SBlaise Bourdin   PetscCheck(snes->usersetbounds, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must set SNESVI bounds before calling SNESVIGetVariableBounds()");
517cf836535SBlaise Bourdin   if (xl) *xl = snes->xl;
518cf836535SBlaise Bourdin   if (xu) *xu = snes->xu;
5193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
520cf836535SBlaise Bourdin }
521cf836535SBlaise Bourdin 
522d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject)
523d71ae5a4SJacob Faibussowitsch {
5248afaa268SBarry Smith   PetscBool flg = PETSC_FALSE;
5252b4ed282SShri Abhyankar 
5262b4ed282SShri Abhyankar   PetscFunctionBegin;
527d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options");
5289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL));
5299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL));
5301baa6e33SBarry Smith   if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL));
531de34d3e9SBarry Smith   flg = PETSC_FALSE;
5329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_vi_monitor_residual", "Monitor residual all non-active variables; using zero for active constraints", "SNESMonitorVIResidual", flg, &flg, NULL));
5331baa6e33SBarry Smith   if (flg) PetscCall(SNESMonitorSet(snes, SNESVIMonitorResidual, PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)), NULL));
534d0609cedSBarry Smith   PetscOptionsHeadEnd();
5353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
536ed70e96aSJungho Lee }
537