xref: /petsc/src/snes/impls/vi/vi.c (revision 9d1809e211fd477413a86b6864cd4fdcfb2f8f55)
19bd66eb0SPeter Brune #include "viimpl.h"  /*I "petscsnes.h" I*/
2d0af7cd3SBarry Smith 
3d0af7cd3SBarry Smith #undef __FUNCT__
42176524cSBarry Smith #define __FUNCT__ "SNESVISetComputeVariableBounds"
52176524cSBarry Smith /*@C
62176524cSBarry Smith    SNESVISetComputeVariableBounds - Sets a function  that is called to compute the variable bounds
72176524cSBarry Smith 
82176524cSBarry Smith    Input parameter
92176524cSBarry Smith +  snes - the SNES context
102176524cSBarry Smith -  compute - computes the bounds
112176524cSBarry Smith 
122bd2b0e6SSatish Balay    Level: advanced
132176524cSBarry Smith 
14aab9d709SJed Brown @*/
1577cdaf69SJed Brown PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
162176524cSBarry Smith {
1761589011SJed Brown   PetscErrorCode   ierr,(*f)(SNES,PetscErrorCode(*)(SNES,Vec,Vec));
182176524cSBarry Smith 
192176524cSBarry Smith   PetscFunctionBegin;
2061589011SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2161589011SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
2261589011SJed Brown   if (!f) {ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);}
2361589011SJed Brown   ierr = PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode(*)(SNES,Vec,Vec)),(snes,compute));CHKERRQ(ierr);
242176524cSBarry Smith   PetscFunctionReturn(0);
252176524cSBarry Smith }
262176524cSBarry Smith 
2761589011SJed Brown EXTERN_C_BEGIN
2861589011SJed Brown #undef __FUNCT__
2961589011SJed Brown #define __FUNCT__ "SNESVISetComputeVariableBounds_VI"
3061589011SJed Brown PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
3161589011SJed Brown {
3261589011SJed Brown   PetscFunctionBegin;
3361589011SJed Brown   snes->ops->computevariablebounds = compute;
3461589011SJed Brown   PetscFunctionReturn(0);
3561589011SJed Brown }
3661589011SJed Brown EXTERN_C_END
372176524cSBarry Smith 
382176524cSBarry Smith #undef __FUNCT__
39ed70e96aSJungho Lee #define __FUNCT__ "SNESVIComputeInactiveSetIS"
40d0af7cd3SBarry Smith /*
41ed70e96aSJungho Lee    SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables
42d0af7cd3SBarry Smith 
43d0af7cd3SBarry Smith    Input parameter
44d0af7cd3SBarry Smith .  snes - the SNES context
45d0af7cd3SBarry Smith .  X    - the snes solution vector
46d0af7cd3SBarry Smith 
47d0af7cd3SBarry Smith    Output parameter
48d0af7cd3SBarry Smith .  ISact - active set index set
49d0af7cd3SBarry Smith 
50d0af7cd3SBarry Smith  */
51ed70e96aSJungho Lee PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS* inact)
52d0af7cd3SBarry Smith {
53d0af7cd3SBarry Smith   PetscErrorCode    ierr;
54d0af7cd3SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
55d0af7cd3SBarry Smith   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
56d0af7cd3SBarry Smith 
57d0af7cd3SBarry Smith   PetscFunctionBegin;
58d0af7cd3SBarry Smith   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
59d0af7cd3SBarry Smith   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
60d0af7cd3SBarry Smith   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
61d0af7cd3SBarry Smith   ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr);
62d0af7cd3SBarry Smith   ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr);
63d0af7cd3SBarry Smith   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
64d0af7cd3SBarry Smith   /* Compute inactive set size */
65d0af7cd3SBarry Smith   for (i=0; i < nlocal;i++) {
66d0af7cd3SBarry Smith     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) nloc_isact++;
67d0af7cd3SBarry Smith   }
68d0af7cd3SBarry Smith 
69d0af7cd3SBarry Smith   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
70d0af7cd3SBarry Smith 
71d0af7cd3SBarry Smith   /* Set inactive set indices */
72d0af7cd3SBarry Smith   for(i=0; i < nlocal; i++) {
73d0af7cd3SBarry Smith     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) idx_act[i1++] = ilow+i;
74d0af7cd3SBarry Smith   }
75d0af7cd3SBarry Smith 
76d0af7cd3SBarry Smith    /* Create inactive set IS */
77d0af7cd3SBarry Smith   ierr = ISCreateGeneral(((PetscObject)upper)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr);
78d0af7cd3SBarry Smith 
79d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
80d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr);
81d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr);
82d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
83d0af7cd3SBarry Smith   PetscFunctionReturn(0);
84d0af7cd3SBarry Smith }
85d0af7cd3SBarry Smith 
863c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/
872b4ed282SShri Abhyankar 
889308c137SBarry Smith #undef __FUNCT__
899308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI"
907087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
919308c137SBarry Smith {
929308c137SBarry Smith   PetscErrorCode    ierr;
93649052a6SBarry Smith   PetscViewer        viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
949308c137SBarry Smith   const PetscScalar  *x,*xl,*xu,*f;
956fd67740SBarry Smith   PetscInt           i,n,act[2] = {0,0},fact[2],N;
966a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
976a9e2d46SJungho Lee   PetscInt           act_bound[2] = {0,0},fact_bound[2];
989308c137SBarry Smith   PetscReal          rnorm,fnorm;
99*9d1809e2SSatish Balay   double             tmp;
1009308c137SBarry Smith 
1019308c137SBarry Smith   PetscFunctionBegin;
1029308c137SBarry Smith   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
1036fd67740SBarry Smith   ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr);
104c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
105c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
1069308c137SBarry Smith   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
1073f731af5SBarry Smith   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
1089308c137SBarry Smith 
1099308c137SBarry Smith   rnorm = 0.0;
1109308c137SBarry Smith   for (i=0; i<n; i++) {
11110f5ae6bSBarry Smith     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
112e400df7aSBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++;
113e400df7aSBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++;
114e400df7aSBarry Smith     else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Can never get here");
1159308c137SBarry Smith   }
1166a9e2d46SJungho Lee 
1176a9e2d46SJungho Lee   for (i=0; i<n; i++) {
1186a9e2d46SJungho Lee     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++;
1196a9e2d46SJungho Lee     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++;
1206a9e2d46SJungho Lee   }
1213f731af5SBarry Smith   ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
122c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
123c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
1249308c137SBarry Smith   ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
125d9822059SBarry Smith   ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
12621a4c9b1SBarry Smith   ierr = MPI_Allreduce(act,fact,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
1276a9e2d46SJungho Lee   ierr = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
128f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1296fd67740SBarry Smith 
130649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
131*9d1809e2SSatish Balay   if(snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds);
132*9d1809e2SSatish Balay   else tmp = 0.0;
133*9d1809e2SSatish Balay   ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES VI Function norm %14.12e Active lower constraints %D/%D upper constraints %D/%D 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);CHKERRQ(ierr);
1346a9e2d46SJungho Lee 
135649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1369308c137SBarry Smith   PetscFunctionReturn(0);
1379308c137SBarry Smith }
1389308c137SBarry Smith 
1392b4ed282SShri Abhyankar /*
1402b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1412b4ed282SShri 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
1422b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1432b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1442b4ed282SShri Abhyankar */
1452b4ed282SShri Abhyankar #undef __FUNCT__
1462b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private"
147ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
1482b4ed282SShri Abhyankar {
1492b4ed282SShri Abhyankar   PetscReal      a1;
1502b4ed282SShri Abhyankar   PetscErrorCode ierr;
151ace3abfcSBarry Smith   PetscBool     hastranspose;
1522b4ed282SShri Abhyankar 
1532b4ed282SShri Abhyankar   PetscFunctionBegin;
1542b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1552b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1562b4ed282SShri Abhyankar   if (hastranspose) {
1572b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1582b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
1592b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
1604839bfe8SBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
1612b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
1622b4ed282SShri Abhyankar   } else {
1632b4ed282SShri Abhyankar     Vec         work;
1642b4ed282SShri Abhyankar     PetscScalar result;
1652b4ed282SShri Abhyankar     PetscReal   wnorm;
1662b4ed282SShri Abhyankar 
1672b4ed282SShri Abhyankar     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
1682b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
1692b4ed282SShri Abhyankar     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
1702b4ed282SShri Abhyankar     ierr = MatMult(A,W,work);CHKERRQ(ierr);
1712b4ed282SShri Abhyankar     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
1726bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
1732b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
1744839bfe8SBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
1752b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1762b4ed282SShri Abhyankar   }
1772b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1782b4ed282SShri Abhyankar }
1792b4ed282SShri Abhyankar 
1802b4ed282SShri Abhyankar /*
1812b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
1822b4ed282SShri Abhyankar */
1832b4ed282SShri Abhyankar #undef __FUNCT__
1842b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private"
1852b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
1862b4ed282SShri Abhyankar {
1872b4ed282SShri Abhyankar   PetscReal      a1,a2;
1882b4ed282SShri Abhyankar   PetscErrorCode ierr;
189ace3abfcSBarry Smith   PetscBool     hastranspose;
1902b4ed282SShri Abhyankar 
1912b4ed282SShri Abhyankar   PetscFunctionBegin;
1922b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1932b4ed282SShri Abhyankar   if (hastranspose) {
1942b4ed282SShri Abhyankar     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
1952b4ed282SShri Abhyankar     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
1962b4ed282SShri Abhyankar 
1972b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1982b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
1992b4ed282SShri Abhyankar     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
2002b4ed282SShri Abhyankar     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
2012b4ed282SShri Abhyankar     if (a1 != 0.0) {
2024839bfe8SBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
2032b4ed282SShri Abhyankar     }
2042b4ed282SShri Abhyankar   }
2052b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2062b4ed282SShri Abhyankar }
2072b4ed282SShri Abhyankar 
2082b4ed282SShri Abhyankar /*
2092b4ed282SShri Abhyankar   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
2102b4ed282SShri Abhyankar 
2112b4ed282SShri Abhyankar   Notes:
2122b4ed282SShri Abhyankar   The convergence criterion currently implemented is
2132b4ed282SShri Abhyankar   merit < abstol
2142b4ed282SShri Abhyankar   merit < rtol*merit_initial
2152b4ed282SShri Abhyankar */
2162b4ed282SShri Abhyankar #undef __FUNCT__
2172b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI"
2187fe79bc4SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
2192b4ed282SShri Abhyankar {
2202b4ed282SShri Abhyankar   PetscErrorCode ierr;
2212b4ed282SShri Abhyankar 
2222b4ed282SShri Abhyankar   PetscFunctionBegin;
2232b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2242b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
2252b4ed282SShri Abhyankar 
2262b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
2272b4ed282SShri Abhyankar 
2282b4ed282SShri Abhyankar   if (!it) {
2292b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
2307fe79bc4SShri Abhyankar     snes->ttol = fnorm*snes->rtol;
2312b4ed282SShri Abhyankar   }
2327fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
2332b4ed282SShri Abhyankar     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
2342b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
2357fe79bc4SShri Abhyankar   } else if (fnorm < snes->abstol) {
2364839bfe8SBarry Smith     ierr = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
2372b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
2382b4ed282SShri Abhyankar   } else if (snes->nfuncs >= snes->max_funcs) {
2392b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
2402b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
2412b4ed282SShri Abhyankar   }
2422b4ed282SShri Abhyankar 
2432b4ed282SShri Abhyankar   if (it && !*reason) {
2447fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2454839bfe8SBarry Smith       ierr = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
2462b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2472b4ed282SShri Abhyankar     }
2482b4ed282SShri Abhyankar   }
2492b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2502b4ed282SShri Abhyankar }
2512b4ed282SShri Abhyankar 
252ee657d29SShri Abhyankar 
253c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
254c1a5e521SShri Abhyankar /*
255c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
256c1a5e521SShri Abhyankar 
257c1a5e521SShri Abhyankar    Input Parameters:
258c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
259c1a5e521SShri Abhyankar 
260c1a5e521SShri Abhyankar    Output Parameters:
261c1a5e521SShri Abhyankar .  X - Bound projected X
262c1a5e521SShri Abhyankar 
263c1a5e521SShri Abhyankar */
264c1a5e521SShri Abhyankar 
265c1a5e521SShri Abhyankar #undef __FUNCT__
266c1a5e521SShri Abhyankar #define __FUNCT__ "SNESVIProjectOntoBounds"
267c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
268c1a5e521SShri Abhyankar {
269c1a5e521SShri Abhyankar   PetscErrorCode    ierr;
270c1a5e521SShri Abhyankar   const PetscScalar *xl,*xu;
271c1a5e521SShri Abhyankar   PetscScalar       *x;
272c1a5e521SShri Abhyankar   PetscInt          i,n;
273c1a5e521SShri Abhyankar 
274c1a5e521SShri Abhyankar   PetscFunctionBegin;
275c1a5e521SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
276c1a5e521SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
277c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
278c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
279c1a5e521SShri Abhyankar 
280c1a5e521SShri Abhyankar   for(i = 0;i<n;i++) {
28110f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
28210f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
283c1a5e521SShri Abhyankar   }
284c1a5e521SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
285c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
286c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
287c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
288c1a5e521SShri Abhyankar }
289c1a5e521SShri Abhyankar 
29069c03620SShri Abhyankar 
29169c03620SShri Abhyankar #undef __FUNCT__
292693ddf92SShri Abhyankar #define __FUNCT__ "SNESVIGetActiveSetIS"
293693ddf92SShri Abhyankar /*
294693ddf92SShri Abhyankar    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
295693ddf92SShri Abhyankar 
296693ddf92SShri Abhyankar    Input parameter
297693ddf92SShri Abhyankar .  snes - the SNES context
298693ddf92SShri Abhyankar .  X    - the snes solution vector
299693ddf92SShri Abhyankar .  F    - the nonlinear function vector
300693ddf92SShri Abhyankar 
301693ddf92SShri Abhyankar    Output parameter
302693ddf92SShri Abhyankar .  ISact - active set index set
303693ddf92SShri Abhyankar  */
304693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS* ISact)
305d950fb63SShri Abhyankar {
306d950fb63SShri Abhyankar   PetscErrorCode   ierr;
307c2fc9fa9SBarry Smith   Vec               Xl=snes->xl,Xu=snes->xu;
308693ddf92SShri Abhyankar   const PetscScalar *x,*f,*xl,*xu;
309693ddf92SShri Abhyankar   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
310d950fb63SShri Abhyankar 
311d950fb63SShri Abhyankar   PetscFunctionBegin;
312d950fb63SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
313d950fb63SShri Abhyankar   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
314fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
315fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
316fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
317fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
318693ddf92SShri Abhyankar   /* Compute active set size */
319d950fb63SShri Abhyankar   for (i=0; i < nlocal;i++) {
32010f5ae6bSBarry Smith     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) nloc_isact++;
321d950fb63SShri Abhyankar   }
322693ddf92SShri Abhyankar 
323d950fb63SShri Abhyankar   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
324d950fb63SShri Abhyankar 
325693ddf92SShri Abhyankar   /* Set active set indices */
326d950fb63SShri Abhyankar   for(i=0; i < nlocal; i++) {
32710f5ae6bSBarry Smith     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) idx_act[i1++] = ilow+i;
328d950fb63SShri Abhyankar   }
329d950fb63SShri Abhyankar 
330693ddf92SShri Abhyankar    /* Create active set IS */
33162298a1eSBarry Smith   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
332d950fb63SShri Abhyankar 
333fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
334fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
335fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
336fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
337d950fb63SShri Abhyankar   PetscFunctionReturn(0);
338d950fb63SShri Abhyankar }
339d950fb63SShri Abhyankar 
340693ddf92SShri Abhyankar #undef __FUNCT__
341693ddf92SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS"
342693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS* ISact,IS* ISinact)
343693ddf92SShri Abhyankar {
344693ddf92SShri Abhyankar   PetscErrorCode     ierr;
345693ddf92SShri Abhyankar 
346693ddf92SShri Abhyankar   PetscFunctionBegin;
347693ddf92SShri Abhyankar   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
348693ddf92SShri Abhyankar   ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr);
349693ddf92SShri Abhyankar   PetscFunctionReturn(0);
350693ddf92SShri Abhyankar }
351693ddf92SShri Abhyankar 
352fe835674SShri Abhyankar #undef __FUNCT__
353fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
35410f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm)
355fe835674SShri Abhyankar {
356fe835674SShri Abhyankar   PetscErrorCode    ierr;
357fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
358fe835674SShri Abhyankar   PetscInt          i,n;
359fe835674SShri Abhyankar   PetscReal         rnorm;
360fe835674SShri Abhyankar 
361fe835674SShri Abhyankar   PetscFunctionBegin;
362fe835674SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
363c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
364c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
365fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
366fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
367fe835674SShri Abhyankar   rnorm = 0.0;
368fe835674SShri Abhyankar   for (i=0; i<n; i++) {
36910f5ae6bSBarry Smith     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
370fe835674SShri Abhyankar   }
371fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
372c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
373c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
374fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
375d9822059SBarry Smith   ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
37662d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
377fe835674SShri Abhyankar   PetscFunctionReturn(0);
378fe835674SShri Abhyankar }
379fe835674SShri Abhyankar 
38008da532bSDmitry Karpeev #undef __FUNCT__
38108da532bSDmitry Karpeev #define __FUNCT__ "SNESVIDMComputeVariableBounds"
38208da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
38308da532bSDmitry Karpeev {
38408da532bSDmitry Karpeev   PetscErrorCode     ierr;
38508da532bSDmitry Karpeev   PetscFunctionBegin;
38608da532bSDmitry Karpeev   ierr = DMComputeVariableBounds(snes->dm, xl, xu); CHKERRQ(ierr);
38708da532bSDmitry Karpeev   PetscFunctionReturn(0);
38808da532bSDmitry Karpeev }
38908da532bSDmitry Karpeev 
3902f63a38cSShri Abhyankar 
3912b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3922b4ed282SShri Abhyankar /*
393c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3942b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3952b4ed282SShri Abhyankar 
3962b4ed282SShri Abhyankar    Input Parameter:
3972b4ed282SShri Abhyankar .  snes - the SNES context
3982b4ed282SShri Abhyankar .  x - the solution vector
3992b4ed282SShri Abhyankar 
4002b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
4012b4ed282SShri Abhyankar 
4022b4ed282SShri Abhyankar    Notes:
4032b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
4042b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
4052b4ed282SShri Abhyankar    the call to SNESSolve().
4062b4ed282SShri Abhyankar  */
4072b4ed282SShri Abhyankar #undef __FUNCT__
4082b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
4092b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
4102b4ed282SShri Abhyankar {
4112b4ed282SShri Abhyankar   PetscErrorCode ierr;
4122b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
4132b4ed282SShri Abhyankar 
4142b4ed282SShri Abhyankar   PetscFunctionBegin;
41558c9b817SLisandro Dalcin 
4169bd66eb0SPeter Brune   ierr = SNESDefaultGetWork(snes,1);CHKERRQ(ierr);
4176cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
4182b4ed282SShri Abhyankar 
41908da532bSDmitry Karpeev   if(!snes->ops->computevariablebounds && snes->dm) {
42008da532bSDmitry Karpeev     snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
42108da532bSDmitry Karpeev   }
422c2fc9fa9SBarry Smith   if (snes->ops->computevariablebounds) {
423c2fc9fa9SBarry Smith     if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
424c2fc9fa9SBarry Smith     if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
425c2fc9fa9SBarry Smith     ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
426c2fc9fa9SBarry Smith   } else if (!snes->xl && !snes->xu) {
4272176524cSBarry Smith     /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
428c2fc9fa9SBarry Smith     ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr);
429c2fc9fa9SBarry Smith     ierr = VecSet(snes->xl,SNES_VI_NINF);CHKERRQ(ierr);
430c2fc9fa9SBarry Smith     ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr);
431c2fc9fa9SBarry Smith     ierr = VecSet(snes->xu,SNES_VI_INF);CHKERRQ(ierr);
4322b4ed282SShri Abhyankar   } else {
4332b4ed282SShri Abhyankar     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
4342b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
435c2fc9fa9SBarry Smith     ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr);
436c2fc9fa9SBarry Smith     ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr);
4372b4ed282SShri 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]))
4382b4ed282SShri 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.");
4392b4ed282SShri Abhyankar   }
440abcba341SShri Abhyankar 
4412b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4422b4ed282SShri Abhyankar }
4432b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
4442176524cSBarry Smith #undef __FUNCT__
4452176524cSBarry Smith #define __FUNCT__ "SNESReset_VI"
4462176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes)
4472176524cSBarry Smith {
4482176524cSBarry Smith   PetscErrorCode ierr;
4492176524cSBarry Smith 
4502176524cSBarry Smith   PetscFunctionBegin;
451c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
452c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
4532176524cSBarry Smith   PetscFunctionReturn(0);
4542176524cSBarry Smith }
4552176524cSBarry Smith 
4562b4ed282SShri Abhyankar /*
4572b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
4582b4ed282SShri Abhyankar    with SNESCreate_VI().
4592b4ed282SShri Abhyankar 
4602b4ed282SShri Abhyankar    Input Parameter:
4612b4ed282SShri Abhyankar .  snes - the SNES context
4622b4ed282SShri Abhyankar 
4632b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4642b4ed282SShri Abhyankar  */
4652b4ed282SShri Abhyankar #undef __FUNCT__
4662b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI"
4672b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
4682b4ed282SShri Abhyankar {
4692b4ed282SShri Abhyankar   PetscErrorCode ierr;
4702b4ed282SShri Abhyankar 
4712b4ed282SShri Abhyankar   PetscFunctionBegin;
4722b4ed282SShri Abhyankar   ierr = PetscFree(snes->data);CHKERRQ(ierr);
4732b4ed282SShri Abhyankar 
4742b4ed282SShri Abhyankar   /* clear composed functions */
4752b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
476c92abb8aSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
4772b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4782b4ed282SShri Abhyankar }
4797fe79bc4SShri Abhyankar 
4805559b345SBarry Smith #undef __FUNCT__
4815559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds"
4825559b345SBarry Smith /*@
4832b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
4842b4ed282SShri Abhyankar 
4852b4ed282SShri Abhyankar    Input Parameters:
4862b4ed282SShri Abhyankar .  snes - the SNES context.
4872b4ed282SShri Abhyankar .  xl   - lower bound.
4882b4ed282SShri Abhyankar .  xu   - upper bound.
4892b4ed282SShri Abhyankar 
4902b4ed282SShri Abhyankar    Notes:
4912b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
49201f0b76bSBarry Smith    SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp().
49384c105d7SBarry Smith 
4942bd2b0e6SSatish Balay    Level: advanced
4952bd2b0e6SSatish Balay 
4965559b345SBarry Smith @*/
49769c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
4982b4ed282SShri Abhyankar {
49961589011SJed Brown   PetscErrorCode   ierr,(*f)(SNES,Vec,Vec);
5002b4ed282SShri Abhyankar 
5012b4ed282SShri Abhyankar   PetscFunctionBegin;
5022b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
5032b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
5042b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
50561589011SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
50661589011SJed Brown   if (!f) {ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);}
50761589011SJed Brown   ierr = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr);
50861589011SJed Brown   PetscFunctionReturn(0);
50961589011SJed Brown }
51061589011SJed Brown 
51161589011SJed Brown EXTERN_C_BEGIN
51261589011SJed Brown #undef __FUNCT__
51361589011SJed Brown #define __FUNCT__ "SNESVISetVariableBounds_VI"
51461589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
51561589011SJed Brown {
51661589011SJed Brown   PetscErrorCode    ierr;
51761589011SJed Brown   const PetscScalar *xxl,*xxu;
51861589011SJed Brown   PetscInt          i,n, cnt = 0;
51961589011SJed Brown 
52061589011SJed Brown   PetscFunctionBegin;
521a63bb30eSJed Brown   ierr = SNESGetFunction(snes,&snes->vec_func,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
522a63bb30eSJed Brown   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
5232b4ed282SShri Abhyankar   if (xl->map->N != snes->vec_func->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xl->map->N,snes->vec_func->map->N);
5242b4ed282SShri Abhyankar   if (xu->map->N != snes->vec_func->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xu->map->N,snes->vec_func->map->N);
525c2fc9fa9SBarry Smith   ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);
5262176524cSBarry Smith   ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
5272176524cSBarry Smith   ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
528c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
529c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
530c2fc9fa9SBarry Smith   snes->xl = xl;
531c2fc9fa9SBarry Smith   snes->xu = xu;
5326fd67740SBarry Smith   ierr = VecGetLocalSize(xl,&n);CHKERRQ(ierr);
5336fd67740SBarry Smith   ierr = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr);
5346fd67740SBarry Smith   ierr = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr);
5356fd67740SBarry Smith   for (i=0; i<n; i++) {
5366fd67740SBarry Smith     cnt += ((xxl[i] != SNES_VI_NINF) || (xxu[i] != SNES_VI_INF));
5376fd67740SBarry Smith   }
538c2fc9fa9SBarry Smith   ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
5396fd67740SBarry Smith   ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr);
5406fd67740SBarry Smith   ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr);
5412b4ed282SShri Abhyankar   PetscFunctionReturn(0);
5422b4ed282SShri Abhyankar }
54361589011SJed Brown EXTERN_C_END
54492c02d66SPeter Brune 
5452b4ed282SShri Abhyankar #undef __FUNCT__
546c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VI"
547c2fc9fa9SBarry Smith PetscErrorCode SNESSetFromOptions_VI(SNES snes)
5482b4ed282SShri Abhyankar {
5492b4ed282SShri Abhyankar   PetscErrorCode  ierr;
550c2fc9fa9SBarry Smith   PetscBool       flg;
5512b4ed282SShri Abhyankar 
5522b4ed282SShri Abhyankar   PetscFunctionBegin;
553c2fc9fa9SBarry Smith   ierr = PetscOptionsHead("SNES VI options");CHKERRQ(ierr);
554c2fc9fa9SBarry Smith   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
555c2fc9fa9SBarry Smith   if (flg) {
556c2fc9fa9SBarry Smith     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
5572b4ed282SShri Abhyankar   }
558c2fc9fa9SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
559ed70e96aSJungho Lee   PetscFunctionReturn(0);
560ed70e96aSJungho Lee }
561