xref: /petsc/src/snes/impls/vi/vi.c (revision 785e854f82a3c614b452fca2cf5ad4f2afe8bdde)
1639ea3faSPeter Brune #include <petsc-private/snesimpl.h>  /*I "petscsnes.h" I*/
21e25c274SJed Brown #include <petscdm.h>
3d0af7cd3SBarry Smith 
4d0af7cd3SBarry Smith #undef __FUNCT__
52176524cSBarry Smith #define __FUNCT__ "SNESVISetComputeVariableBounds"
62176524cSBarry Smith /*@C
72176524cSBarry Smith    SNESVISetComputeVariableBounds - Sets a function  that is called to compute the variable bounds
82176524cSBarry Smith 
92176524cSBarry Smith    Input parameter
102176524cSBarry Smith +  snes - the SNES context
112176524cSBarry Smith -  compute - computes the bounds
122176524cSBarry Smith 
132bd2b0e6SSatish Balay    Level: advanced
142176524cSBarry Smith 
15c1c3a0ecSBarry Smith .seealso:   SNESVISetVariableBounds()
16c1c3a0ecSBarry Smith 
17aab9d709SJed Brown @*/
1877cdaf69SJed Brown PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
192176524cSBarry Smith {
2061589011SJed Brown   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec));
212176524cSBarry Smith 
222176524cSBarry Smith   PetscFunctionBegin;
2361589011SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
240005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",&f);CHKERRQ(ierr);
25f450aa47SBarry Smith   if (!f) {ierr = SNESSetType(snes,SNESVINEWTONRSLS);CHKERRQ(ierr);}
2661589011SJed Brown   ierr = PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec)),(snes,compute));CHKERRQ(ierr);
272176524cSBarry Smith   PetscFunctionReturn(0);
282176524cSBarry Smith }
292176524cSBarry Smith 
3061589011SJed Brown #undef __FUNCT__
3161589011SJed Brown #define __FUNCT__ "SNESVISetComputeVariableBounds_VI"
3261589011SJed Brown PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
3361589011SJed Brown {
3461589011SJed Brown   PetscFunctionBegin;
3561589011SJed Brown   snes->ops->computevariablebounds = compute;
3661589011SJed Brown   PetscFunctionReturn(0);
3761589011SJed Brown }
382176524cSBarry Smith 
392176524cSBarry Smith #undef __FUNCT__
40ed70e96aSJungho Lee #define __FUNCT__ "SNESVIComputeInactiveSetIS"
41d0af7cd3SBarry Smith /*
42ed70e96aSJungho Lee    SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables
43d0af7cd3SBarry Smith 
44d0af7cd3SBarry Smith    Input parameter
45d0af7cd3SBarry Smith .  snes - the SNES context
46d0af7cd3SBarry Smith .  X    - the snes solution vector
47d0af7cd3SBarry Smith 
48d0af7cd3SBarry Smith    Output parameter
49d0af7cd3SBarry Smith .  ISact - active set index set
50d0af7cd3SBarry Smith 
51d0af7cd3SBarry Smith  */
52ed70e96aSJungho Lee PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS *inact)
53d0af7cd3SBarry Smith {
54d0af7cd3SBarry Smith   PetscErrorCode    ierr;
55d0af7cd3SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
56d0af7cd3SBarry Smith   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
57d0af7cd3SBarry Smith 
58d0af7cd3SBarry Smith   PetscFunctionBegin;
59d0af7cd3SBarry Smith   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
60d0af7cd3SBarry Smith   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
61d0af7cd3SBarry Smith   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
62d0af7cd3SBarry Smith   ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr);
63d0af7cd3SBarry Smith   ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr);
64d0af7cd3SBarry Smith   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
65d0af7cd3SBarry Smith   /* Compute inactive set size */
66d0af7cd3SBarry Smith   for (i=0; i < nlocal; i++) {
67d0af7cd3SBarry 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++;
68d0af7cd3SBarry Smith   }
69d0af7cd3SBarry Smith 
70*785e854fSJed Brown   ierr = PetscMalloc1(nloc_isact,&idx_act);CHKERRQ(ierr);
71d0af7cd3SBarry Smith 
72d0af7cd3SBarry Smith   /* Set inactive set indices */
73d0af7cd3SBarry Smith   for (i=0; i < nlocal; i++) {
74d0af7cd3SBarry 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;
75d0af7cd3SBarry Smith   }
76d0af7cd3SBarry Smith 
77d0af7cd3SBarry Smith   /* Create inactive set IS */
78ce94432eSBarry Smith   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)upper),nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr);
79d0af7cd3SBarry Smith 
80d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
81d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr);
82d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr);
83d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
84d0af7cd3SBarry Smith   PetscFunctionReturn(0);
85d0af7cd3SBarry Smith }
86d0af7cd3SBarry Smith 
873c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/
882b4ed282SShri Abhyankar 
899308c137SBarry Smith #undef __FUNCT__
909308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI"
917087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
929308c137SBarry Smith {
939308c137SBarry Smith   PetscErrorCode    ierr;
94ce94432eSBarry Smith   PetscViewer       viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
959308c137SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
966fd67740SBarry Smith   PetscInt          i,n,act[2] = {0,0},fact[2],N;
976a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
986a9e2d46SJungho Lee   PetscInt  act_bound[2] = {0,0},fact_bound[2];
999308c137SBarry Smith   PetscReal rnorm,fnorm;
1009d1809e2SSatish Balay   double    tmp;
1019308c137SBarry Smith 
1029308c137SBarry Smith   PetscFunctionBegin;
1039308c137SBarry Smith   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
1046fd67740SBarry Smith   ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr);
105c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
106c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
1079308c137SBarry Smith   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
1083f731af5SBarry Smith   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
1099308c137SBarry Smith 
1109308c137SBarry Smith   rnorm = 0.0;
1119308c137SBarry Smith   for (i=0; i<n; i++) {
11210f5ae6bSBarry 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]);
113e400df7aSBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++;
114e400df7aSBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++;
115ce94432eSBarry Smith     else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Can never get here");
1169308c137SBarry Smith   }
1176a9e2d46SJungho Lee 
1186a9e2d46SJungho Lee   for (i=0; i<n; i++) {
1196a9e2d46SJungho Lee     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++;
1206a9e2d46SJungho Lee     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++;
1216a9e2d46SJungho Lee   }
1223f731af5SBarry Smith   ierr  = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
123c2fc9fa9SBarry Smith   ierr  = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
124c2fc9fa9SBarry Smith   ierr  = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
1259308c137SBarry Smith   ierr  = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
126ce94432eSBarry Smith   ierr  = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
127ce94432eSBarry Smith   ierr  = MPI_Allreduce(act,fact,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
128ce94432eSBarry Smith   ierr  = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
129f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1306fd67740SBarry Smith 
131649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1329d1809e2SSatish Balay   if (snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds);
1339d1809e2SSatish Balay   else tmp = 0.0;
1349d1809e2SSatish 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);
1356a9e2d46SJungho Lee 
136649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1379308c137SBarry Smith   PetscFunctionReturn(0);
1389308c137SBarry Smith }
1399308c137SBarry Smith 
1402b4ed282SShri Abhyankar /*
1412b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1422b4ed282SShri 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
1432b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1442b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1452b4ed282SShri Abhyankar */
1462b4ed282SShri Abhyankar #undef __FUNCT__
1472b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private"
148ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
1492b4ed282SShri Abhyankar {
1502b4ed282SShri Abhyankar   PetscReal      a1;
1512b4ed282SShri Abhyankar   PetscErrorCode ierr;
152ace3abfcSBarry Smith   PetscBool      hastranspose;
1532b4ed282SShri Abhyankar 
1542b4ed282SShri Abhyankar   PetscFunctionBegin;
1552b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1562b4ed282SShri Abhyankar   ierr   = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1572b4ed282SShri Abhyankar   if (hastranspose) {
1582b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1592b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
1602b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
1614839bfe8SBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
1622b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
1632b4ed282SShri Abhyankar   } else {
1642b4ed282SShri Abhyankar     Vec         work;
1652b4ed282SShri Abhyankar     PetscScalar result;
1662b4ed282SShri Abhyankar     PetscReal   wnorm;
1672b4ed282SShri Abhyankar 
1680298fd71SBarry Smith     ierr = VecSetRandom(W,NULL);CHKERRQ(ierr);
1692b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
1702b4ed282SShri Abhyankar     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
1712b4ed282SShri Abhyankar     ierr = MatMult(A,W,work);CHKERRQ(ierr);
1722b4ed282SShri Abhyankar     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
1736bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
1742b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
1754839bfe8SBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
1762b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1772b4ed282SShri Abhyankar   }
1782b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1792b4ed282SShri Abhyankar }
1802b4ed282SShri Abhyankar 
1812b4ed282SShri Abhyankar /*
1822b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
1832b4ed282SShri Abhyankar */
1842b4ed282SShri Abhyankar #undef __FUNCT__
1852b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private"
1862b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
1872b4ed282SShri Abhyankar {
1882b4ed282SShri Abhyankar   PetscReal      a1,a2;
1892b4ed282SShri Abhyankar   PetscErrorCode ierr;
190ace3abfcSBarry Smith   PetscBool      hastranspose;
1912b4ed282SShri Abhyankar 
1922b4ed282SShri Abhyankar   PetscFunctionBegin;
1932b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1942b4ed282SShri Abhyankar   if (hastranspose) {
1952b4ed282SShri Abhyankar     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
1962b4ed282SShri Abhyankar     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
1972b4ed282SShri Abhyankar 
1982b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1992b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
2002b4ed282SShri Abhyankar     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
2012b4ed282SShri Abhyankar     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
2022b4ed282SShri Abhyankar     if (a1 != 0.0) {
2034839bfe8SBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
2042b4ed282SShri Abhyankar     }
2052b4ed282SShri Abhyankar   }
2062b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2072b4ed282SShri Abhyankar }
2082b4ed282SShri Abhyankar 
2092b4ed282SShri Abhyankar /*
2108d359177SBarry Smith   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
2112b4ed282SShri Abhyankar 
2122b4ed282SShri Abhyankar   Notes:
2132b4ed282SShri Abhyankar   The convergence criterion currently implemented is
2142b4ed282SShri Abhyankar   merit < abstol
2152b4ed282SShri Abhyankar   merit < rtol*merit_initial
2162b4ed282SShri Abhyankar */
2172b4ed282SShri Abhyankar #undef __FUNCT__
2188d359177SBarry Smith #define __FUNCT__ "SNESConvergedDefault_VI"
2198d359177SBarry Smith PetscErrorCode SNESConvergedDefault_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
2202b4ed282SShri Abhyankar {
2212b4ed282SShri Abhyankar   PetscErrorCode ierr;
2222b4ed282SShri Abhyankar 
2232b4ed282SShri Abhyankar   PetscFunctionBegin;
2242b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2252b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
2262b4ed282SShri Abhyankar 
2272b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
2282b4ed282SShri Abhyankar 
2292b4ed282SShri Abhyankar   if (!it) {
2302b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
2317fe79bc4SShri Abhyankar     snes->ttol = fnorm*snes->rtol;
2322b4ed282SShri Abhyankar   }
2337fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
2342b4ed282SShri Abhyankar     ierr    = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
2352b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
2367fe79bc4SShri Abhyankar   } else if (fnorm < snes->abstol) {
2374839bfe8SBarry Smith     ierr    = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
2382b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
2392b4ed282SShri Abhyankar   } else if (snes->nfuncs >= snes->max_funcs) {
2402b4ed282SShri Abhyankar     ierr    = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
2412b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
2422b4ed282SShri Abhyankar   }
2432b4ed282SShri Abhyankar 
2442b4ed282SShri Abhyankar   if (it && !*reason) {
2457fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2464839bfe8SBarry Smith       ierr    = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
2472b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2482b4ed282SShri Abhyankar     }
2492b4ed282SShri Abhyankar   }
2502b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2512b4ed282SShri Abhyankar }
2522b4ed282SShri Abhyankar 
253ee657d29SShri Abhyankar 
254c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
255c1a5e521SShri Abhyankar /*
256c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
257c1a5e521SShri Abhyankar 
258c1a5e521SShri Abhyankar    Input Parameters:
259c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
260c1a5e521SShri Abhyankar 
261c1a5e521SShri Abhyankar    Output Parameters:
262c1a5e521SShri Abhyankar .  X - Bound projected X
263c1a5e521SShri Abhyankar 
264c1a5e521SShri Abhyankar */
265c1a5e521SShri Abhyankar 
266c1a5e521SShri Abhyankar #undef __FUNCT__
267c1a5e521SShri Abhyankar #define __FUNCT__ "SNESVIProjectOntoBounds"
268c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
269c1a5e521SShri Abhyankar {
270c1a5e521SShri Abhyankar   PetscErrorCode    ierr;
271c1a5e521SShri Abhyankar   const PetscScalar *xl,*xu;
272c1a5e521SShri Abhyankar   PetscScalar       *x;
273c1a5e521SShri Abhyankar   PetscInt          i,n;
274c1a5e521SShri Abhyankar 
275c1a5e521SShri Abhyankar   PetscFunctionBegin;
276c1a5e521SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
277c1a5e521SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
278c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
279c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
280c1a5e521SShri Abhyankar 
281c1a5e521SShri Abhyankar   for (i = 0; i<n; i++) {
28210f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
28310f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
284c1a5e521SShri Abhyankar   }
285c1a5e521SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
286c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
287c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
288c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
289c1a5e521SShri Abhyankar }
290c1a5e521SShri Abhyankar 
29169c03620SShri Abhyankar 
29269c03620SShri Abhyankar #undef __FUNCT__
293693ddf92SShri Abhyankar #define __FUNCT__ "SNESVIGetActiveSetIS"
294693ddf92SShri Abhyankar /*
295693ddf92SShri Abhyankar    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
296693ddf92SShri Abhyankar 
297693ddf92SShri Abhyankar    Input parameter
298693ddf92SShri Abhyankar .  snes - the SNES context
299693ddf92SShri Abhyankar .  X    - the snes solution vector
300693ddf92SShri Abhyankar .  F    - the nonlinear function vector
301693ddf92SShri Abhyankar 
302693ddf92SShri Abhyankar    Output parameter
303693ddf92SShri Abhyankar .  ISact - active set index set
304693ddf92SShri Abhyankar  */
305693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS *ISact)
306d950fb63SShri Abhyankar {
307d950fb63SShri Abhyankar   PetscErrorCode    ierr;
308c2fc9fa9SBarry Smith   Vec               Xl=snes->xl,Xu=snes->xu;
309693ddf92SShri Abhyankar   const PetscScalar *x,*f,*xl,*xu;
310693ddf92SShri Abhyankar   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
311d950fb63SShri Abhyankar 
312d950fb63SShri Abhyankar   PetscFunctionBegin;
313d950fb63SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
314d950fb63SShri Abhyankar   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
315fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
316fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
317fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
318fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
319693ddf92SShri Abhyankar   /* Compute active set size */
320d950fb63SShri Abhyankar   for (i=0; i < nlocal;i++) {
32110f5ae6bSBarry 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++;
322d950fb63SShri Abhyankar   }
323693ddf92SShri Abhyankar 
324*785e854fSJed Brown   ierr = PetscMalloc1(nloc_isact,&idx_act);CHKERRQ(ierr);
325d950fb63SShri Abhyankar 
326693ddf92SShri Abhyankar   /* Set active set indices */
327d950fb63SShri Abhyankar   for (i=0; i < nlocal; i++) {
32810f5ae6bSBarry 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;
329d950fb63SShri Abhyankar   }
330d950fb63SShri Abhyankar 
331693ddf92SShri Abhyankar   /* Create active set IS */
332ce94432eSBarry Smith   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes),nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
333d950fb63SShri Abhyankar 
334fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
335fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
336fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
337fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
338d950fb63SShri Abhyankar   PetscFunctionReturn(0);
339d950fb63SShri Abhyankar }
340d950fb63SShri Abhyankar 
341693ddf92SShri Abhyankar #undef __FUNCT__
342693ddf92SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS"
343693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact)
344693ddf92SShri Abhyankar {
345693ddf92SShri Abhyankar   PetscErrorCode ierr;
346077aedafSJed Brown   PetscInt       rstart,rend;
347693ddf92SShri Abhyankar 
348693ddf92SShri Abhyankar   PetscFunctionBegin;
349693ddf92SShri Abhyankar   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
350077aedafSJed Brown   ierr = VecGetOwnershipRange(X,&rstart,&rend);CHKERRQ(ierr);
351077aedafSJed Brown   ierr = ISComplement(*ISact,rstart,rend,ISinact);CHKERRQ(ierr);
352693ddf92SShri Abhyankar   PetscFunctionReturn(0);
353693ddf92SShri Abhyankar }
354693ddf92SShri Abhyankar 
355fe835674SShri Abhyankar #undef __FUNCT__
356fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
35710f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm)
358fe835674SShri Abhyankar {
359fe835674SShri Abhyankar   PetscErrorCode    ierr;
360fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
361fe835674SShri Abhyankar   PetscInt          i,n;
362fe835674SShri Abhyankar   PetscReal         rnorm;
363fe835674SShri Abhyankar 
364fe835674SShri Abhyankar   PetscFunctionBegin;
365fe835674SShri Abhyankar   ierr  = VecGetLocalSize(X,&n);CHKERRQ(ierr);
366c2fc9fa9SBarry Smith   ierr  = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
367c2fc9fa9SBarry Smith   ierr  = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
368fe835674SShri Abhyankar   ierr  = VecGetArrayRead(X,&x);CHKERRQ(ierr);
369fe835674SShri Abhyankar   ierr  = VecGetArrayRead(F,&f);CHKERRQ(ierr);
370fe835674SShri Abhyankar   rnorm = 0.0;
371fe835674SShri Abhyankar   for (i=0; i<n; i++) {
37210f5ae6bSBarry 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]);
3738f918527SKarl Rupp   }
374fe835674SShri Abhyankar   ierr   = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
375c2fc9fa9SBarry Smith   ierr   = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
376c2fc9fa9SBarry Smith   ierr   = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
377fe835674SShri Abhyankar   ierr   = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
378ce94432eSBarry Smith   ierr   = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
37962d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
380fe835674SShri Abhyankar   PetscFunctionReturn(0);
381fe835674SShri Abhyankar }
382fe835674SShri Abhyankar 
38308da532bSDmitry Karpeev #undef __FUNCT__
38408da532bSDmitry Karpeev #define __FUNCT__ "SNESVIDMComputeVariableBounds"
38508da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
38608da532bSDmitry Karpeev {
38708da532bSDmitry Karpeev   PetscErrorCode ierr;
3886e111a19SKarl Rupp 
38908da532bSDmitry Karpeev   PetscFunctionBegin;
39008da532bSDmitry Karpeev   ierr = DMComputeVariableBounds(snes->dm, xl, xu);CHKERRQ(ierr);
39108da532bSDmitry Karpeev   PetscFunctionReturn(0);
39208da532bSDmitry Karpeev }
39308da532bSDmitry Karpeev 
3942f63a38cSShri Abhyankar 
3952b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3962b4ed282SShri Abhyankar /*
397c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3982b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3992b4ed282SShri Abhyankar 
4002b4ed282SShri Abhyankar    Input Parameter:
4012b4ed282SShri Abhyankar .  snes - the SNES context
4022b4ed282SShri Abhyankar .  x - the solution vector
4032b4ed282SShri Abhyankar 
4042b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
4052b4ed282SShri Abhyankar 
4062b4ed282SShri Abhyankar    Notes:
4072b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
4082b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
4092b4ed282SShri Abhyankar    the call to SNESSolve().
4102b4ed282SShri Abhyankar  */
4112b4ed282SShri Abhyankar #undef __FUNCT__
4122b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
4132b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
4142b4ed282SShri Abhyankar {
4152b4ed282SShri Abhyankar   PetscErrorCode ierr;
4162b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
4172b4ed282SShri Abhyankar 
4182b4ed282SShri Abhyankar   PetscFunctionBegin;
419fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,1);CHKERRQ(ierr);
4206cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
4212b4ed282SShri Abhyankar 
42208da532bSDmitry Karpeev   if (!snes->ops->computevariablebounds && snes->dm) {
423a201590fSDmitry Karpeev     PetscBool flag;
424a201590fSDmitry Karpeev     ierr = DMHasVariableBounds(snes->dm, &flag);CHKERRQ(ierr);
4251aa26658SKarl Rupp 
42608da532bSDmitry Karpeev     snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
42708da532bSDmitry Karpeev   }
428a201590fSDmitry Karpeev   if (!snes->usersetbounds) {
429c2fc9fa9SBarry Smith     if (snes->ops->computevariablebounds) {
430c2fc9fa9SBarry Smith       if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
431c2fc9fa9SBarry Smith       if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
432c2fc9fa9SBarry Smith       ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
4331aa26658SKarl Rupp     } else if (!snes->xl && !snes->xu) {
4342176524cSBarry Smith       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
435c2fc9fa9SBarry Smith       ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr);
436c2fc9fa9SBarry Smith       ierr = VecSet(snes->xl,SNES_VI_NINF);CHKERRQ(ierr);
437c2fc9fa9SBarry Smith       ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr);
438c2fc9fa9SBarry Smith       ierr = VecSet(snes->xu,SNES_VI_INF);CHKERRQ(ierr);
4392b4ed282SShri Abhyankar     } else {
4402b4ed282SShri Abhyankar       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
4412b4ed282SShri Abhyankar       ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
442c2fc9fa9SBarry Smith       ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr);
443c2fc9fa9SBarry Smith       ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr);
4442b4ed282SShri 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]))
4452b4ed282SShri 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.");
4462b4ed282SShri Abhyankar     }
447a201590fSDmitry Karpeev   }
4482b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4492b4ed282SShri Abhyankar }
4502b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
4512176524cSBarry Smith #undef __FUNCT__
4522176524cSBarry Smith #define __FUNCT__ "SNESReset_VI"
4532176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes)
4542176524cSBarry Smith {
4552176524cSBarry Smith   PetscErrorCode ierr;
4562176524cSBarry Smith 
4572176524cSBarry Smith   PetscFunctionBegin;
458c2fc9fa9SBarry Smith   ierr                = VecDestroy(&snes->xl);CHKERRQ(ierr);
459c2fc9fa9SBarry Smith   ierr                = VecDestroy(&snes->xu);CHKERRQ(ierr);
4602d6615e8SDmitry Karpeev   snes->usersetbounds = PETSC_FALSE;
4612176524cSBarry Smith   PetscFunctionReturn(0);
4622176524cSBarry Smith }
4632176524cSBarry Smith 
4642b4ed282SShri Abhyankar /*
4652b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
4662b4ed282SShri Abhyankar    with SNESCreate_VI().
4672b4ed282SShri Abhyankar 
4682b4ed282SShri Abhyankar    Input Parameter:
4692b4ed282SShri Abhyankar .  snes - the SNES context
4702b4ed282SShri Abhyankar 
4712b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4722b4ed282SShri Abhyankar  */
4732b4ed282SShri Abhyankar #undef __FUNCT__
4742b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI"
4752b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
4762b4ed282SShri Abhyankar {
4772b4ed282SShri Abhyankar   PetscErrorCode ierr;
4782b4ed282SShri Abhyankar 
4792b4ed282SShri Abhyankar   PetscFunctionBegin;
4802b4ed282SShri Abhyankar   ierr = PetscFree(snes->data);CHKERRQ(ierr);
4812b4ed282SShri Abhyankar 
4822b4ed282SShri Abhyankar   /* clear composed functions */
483bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSet_C",NULL);CHKERRQ(ierr);
484bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSetMonitor_C",NULL);CHKERRQ(ierr);
4852b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4862b4ed282SShri Abhyankar }
4877fe79bc4SShri Abhyankar 
4885559b345SBarry Smith #undef __FUNCT__
4895559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds"
4905559b345SBarry Smith /*@
4912b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
4922b4ed282SShri Abhyankar 
4932b4ed282SShri Abhyankar    Input Parameters:
4942b4ed282SShri Abhyankar .  snes - the SNES context.
4952b4ed282SShri Abhyankar .  xl   - lower bound.
4962b4ed282SShri Abhyankar .  xu   - upper bound.
4972b4ed282SShri Abhyankar 
4982b4ed282SShri Abhyankar    Notes:
4992b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
50001f0b76bSBarry Smith    SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp().
50184c105d7SBarry Smith 
5022bd2b0e6SSatish Balay    Level: advanced
5032bd2b0e6SSatish Balay 
5045559b345SBarry Smith @*/
50569c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
5062b4ed282SShri Abhyankar {
50761589011SJed Brown   PetscErrorCode ierr,(*f)(SNES,Vec,Vec);
5082b4ed282SShri Abhyankar 
5092b4ed282SShri Abhyankar   PetscFunctionBegin;
5102b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
5112b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
5122b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
5130005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",&f);CHKERRQ(ierr);
514f450aa47SBarry Smith   if (!f) {ierr = SNESSetType(snes,SNESVINEWTONRSLS);CHKERRQ(ierr);}
51561589011SJed Brown   ierr                = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr);
516a201590fSDmitry Karpeev   snes->usersetbounds = PETSC_TRUE;
51761589011SJed Brown   PetscFunctionReturn(0);
51861589011SJed Brown }
51961589011SJed Brown 
52061589011SJed Brown #undef __FUNCT__
52161589011SJed Brown #define __FUNCT__ "SNESVISetVariableBounds_VI"
52261589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
52361589011SJed Brown {
52461589011SJed Brown   PetscErrorCode    ierr;
52561589011SJed Brown   const PetscScalar *xxl,*xxu;
52661589011SJed Brown   PetscInt          i,n, cnt = 0;
52761589011SJed Brown 
52861589011SJed Brown   PetscFunctionBegin;
5290298fd71SBarry Smith   ierr = SNESGetFunction(snes,&snes->vec_func,NULL,NULL);CHKERRQ(ierr);
530a63bb30eSJed Brown   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
531077aedafSJed Brown   {
532077aedafSJed Brown     PetscInt xlN,xuN,N;
533077aedafSJed Brown     ierr = VecGetSize(xl,&xlN);CHKERRQ(ierr);
534077aedafSJed Brown     ierr = VecGetSize(xu,&xuN);CHKERRQ(ierr);
535077aedafSJed Brown     ierr = VecGetSize(snes->vec_func,&N);CHKERRQ(ierr);
536077aedafSJed Brown     if (xlN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xlN,N);
537077aedafSJed Brown     if (xuN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xuN,N);
538077aedafSJed Brown   }
539f450aa47SBarry Smith   ierr     = SNESSetType(snes,SNESVINEWTONRSLS);CHKERRQ(ierr);
5402176524cSBarry Smith   ierr     = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
5412176524cSBarry Smith   ierr     = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
542c2fc9fa9SBarry Smith   ierr     = VecDestroy(&snes->xl);CHKERRQ(ierr);
543c2fc9fa9SBarry Smith   ierr     = VecDestroy(&snes->xu);CHKERRQ(ierr);
544c2fc9fa9SBarry Smith   snes->xl = xl;
545c2fc9fa9SBarry Smith   snes->xu = xu;
5466fd67740SBarry Smith   ierr     = VecGetLocalSize(xl,&n);CHKERRQ(ierr);
5476fd67740SBarry Smith   ierr     = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr);
5486fd67740SBarry Smith   ierr     = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr);
5491aa26658SKarl Rupp   for (i=0; i<n; i++) cnt += ((xxl[i] != SNES_VI_NINF) || (xxu[i] != SNES_VI_INF));
5501aa26658SKarl Rupp 
551ce94432eSBarry Smith   ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
5526fd67740SBarry Smith   ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr);
5536fd67740SBarry Smith   ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr);
5542b4ed282SShri Abhyankar   PetscFunctionReturn(0);
5552b4ed282SShri Abhyankar }
55692c02d66SPeter Brune 
5572b4ed282SShri Abhyankar #undef __FUNCT__
558c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VI"
559c2fc9fa9SBarry Smith PetscErrorCode SNESSetFromOptions_VI(SNES snes)
5602b4ed282SShri Abhyankar {
5612b4ed282SShri Abhyankar   PetscErrorCode ierr;
562c2fc9fa9SBarry Smith   PetscBool      flg;
563639ea3faSPeter Brune   SNESLineSearch linesearch;
5642b4ed282SShri Abhyankar 
5652b4ed282SShri Abhyankar   PetscFunctionBegin;
566c2fc9fa9SBarry Smith   ierr = PetscOptionsHead("SNES VI options");CHKERRQ(ierr);
567c2fc9fa9SBarry Smith   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
568c2fc9fa9SBarry Smith   if (flg) {
569c2fc9fa9SBarry Smith     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
5702b4ed282SShri Abhyankar   }
571639ea3faSPeter Brune   if (!snes->linesearch) {
5727601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
573639ea3faSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr);
574639ea3faSPeter Brune     ierr = SNESLineSearchBTSetAlpha(linesearch, 0.0);CHKERRQ(ierr);
575639ea3faSPeter Brune   }
576c2fc9fa9SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
577ed70e96aSJungho Lee   PetscFunctionReturn(0);
578ed70e96aSJungho Lee }
579