xref: /petsc/src/snes/impls/vi/vi.c (revision 1e25c2744a56dc779d05d2a7975a7f2272f1aa7b)
1639ea3faSPeter Brune #include <petsc-private/snesimpl.h>  /*I "petscsnes.h" I*/
2*1e25c274SJed 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);
2461589011SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",(PetscVoidStarFunction)&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 EXTERN_C_BEGIN
3161589011SJed Brown #undef __FUNCT__
3261589011SJed Brown #define __FUNCT__ "SNESVISetComputeVariableBounds_VI"
3361589011SJed Brown PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
3461589011SJed Brown {
3561589011SJed Brown   PetscFunctionBegin;
3661589011SJed Brown   snes->ops->computevariablebounds = compute;
3761589011SJed Brown   PetscFunctionReturn(0);
3861589011SJed Brown }
3961589011SJed Brown EXTERN_C_END
402176524cSBarry Smith 
412176524cSBarry Smith #undef __FUNCT__
42ed70e96aSJungho Lee #define __FUNCT__ "SNESVIComputeInactiveSetIS"
43d0af7cd3SBarry Smith /*
44ed70e96aSJungho Lee    SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables
45d0af7cd3SBarry Smith 
46d0af7cd3SBarry Smith    Input parameter
47d0af7cd3SBarry Smith .  snes - the SNES context
48d0af7cd3SBarry Smith .  X    - the snes solution vector
49d0af7cd3SBarry Smith 
50d0af7cd3SBarry Smith    Output parameter
51d0af7cd3SBarry Smith .  ISact - active set index set
52d0af7cd3SBarry Smith 
53d0af7cd3SBarry Smith  */
54ed70e96aSJungho Lee PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS *inact)
55d0af7cd3SBarry Smith {
56d0af7cd3SBarry Smith   PetscErrorCode    ierr;
57d0af7cd3SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
58d0af7cd3SBarry Smith   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
59d0af7cd3SBarry Smith 
60d0af7cd3SBarry Smith   PetscFunctionBegin;
61d0af7cd3SBarry Smith   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
62d0af7cd3SBarry Smith   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
63d0af7cd3SBarry Smith   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
64d0af7cd3SBarry Smith   ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr);
65d0af7cd3SBarry Smith   ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr);
66d0af7cd3SBarry Smith   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
67d0af7cd3SBarry Smith   /* Compute inactive set size */
68d0af7cd3SBarry Smith   for (i=0; i < nlocal; i++) {
69d0af7cd3SBarry 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++;
70d0af7cd3SBarry Smith   }
71d0af7cd3SBarry Smith 
72d0af7cd3SBarry Smith   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
73d0af7cd3SBarry Smith 
74d0af7cd3SBarry Smith   /* Set inactive set indices */
75d0af7cd3SBarry Smith   for (i=0; i < nlocal; i++) {
76d0af7cd3SBarry 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;
77d0af7cd3SBarry Smith   }
78d0af7cd3SBarry Smith 
79d0af7cd3SBarry Smith   /* Create inactive set IS */
80ce94432eSBarry Smith   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)upper),nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr);
81d0af7cd3SBarry Smith 
82d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
83d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr);
84d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr);
85d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
86d0af7cd3SBarry Smith   PetscFunctionReturn(0);
87d0af7cd3SBarry Smith }
88d0af7cd3SBarry Smith 
893c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/
902b4ed282SShri Abhyankar 
919308c137SBarry Smith #undef __FUNCT__
929308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI"
937087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
949308c137SBarry Smith {
959308c137SBarry Smith   PetscErrorCode    ierr;
96ce94432eSBarry Smith   PetscViewer       viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
979308c137SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
986fd67740SBarry Smith   PetscInt          i,n,act[2] = {0,0},fact[2],N;
996a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
1006a9e2d46SJungho Lee   PetscInt  act_bound[2] = {0,0},fact_bound[2];
1019308c137SBarry Smith   PetscReal rnorm,fnorm;
1029d1809e2SSatish Balay   double    tmp;
1039308c137SBarry Smith 
1049308c137SBarry Smith   PetscFunctionBegin;
1059308c137SBarry Smith   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
1066fd67740SBarry Smith   ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr);
107c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
108c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
1099308c137SBarry Smith   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
1103f731af5SBarry Smith   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
1119308c137SBarry Smith 
1129308c137SBarry Smith   rnorm = 0.0;
1139308c137SBarry Smith   for (i=0; i<n; i++) {
11410f5ae6bSBarry 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]);
115e400df7aSBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++;
116e400df7aSBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++;
117ce94432eSBarry Smith     else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Can never get here");
1189308c137SBarry Smith   }
1196a9e2d46SJungho Lee 
1206a9e2d46SJungho Lee   for (i=0; i<n; i++) {
1216a9e2d46SJungho Lee     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++;
1226a9e2d46SJungho Lee     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++;
1236a9e2d46SJungho Lee   }
1243f731af5SBarry Smith   ierr  = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
125c2fc9fa9SBarry Smith   ierr  = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
126c2fc9fa9SBarry Smith   ierr  = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
1279308c137SBarry Smith   ierr  = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
128ce94432eSBarry Smith   ierr  = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
129ce94432eSBarry Smith   ierr  = MPI_Allreduce(act,fact,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
130ce94432eSBarry Smith   ierr  = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
131f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1326fd67740SBarry Smith 
133649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1349d1809e2SSatish Balay   if (snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds);
1359d1809e2SSatish Balay   else tmp = 0.0;
1369d1809e2SSatish 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);
1376a9e2d46SJungho Lee 
138649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1399308c137SBarry Smith   PetscFunctionReturn(0);
1409308c137SBarry Smith }
1419308c137SBarry Smith 
1422b4ed282SShri Abhyankar /*
1432b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1442b4ed282SShri 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
1452b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1462b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1472b4ed282SShri Abhyankar */
1482b4ed282SShri Abhyankar #undef __FUNCT__
1492b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private"
150ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
1512b4ed282SShri Abhyankar {
1522b4ed282SShri Abhyankar   PetscReal      a1;
1532b4ed282SShri Abhyankar   PetscErrorCode ierr;
154ace3abfcSBarry Smith   PetscBool      hastranspose;
1552b4ed282SShri Abhyankar 
1562b4ed282SShri Abhyankar   PetscFunctionBegin;
1572b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1582b4ed282SShri Abhyankar   ierr   = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1592b4ed282SShri Abhyankar   if (hastranspose) {
1602b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1612b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
1622b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
1634839bfe8SBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
1642b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
1652b4ed282SShri Abhyankar   } else {
1662b4ed282SShri Abhyankar     Vec         work;
1672b4ed282SShri Abhyankar     PetscScalar result;
1682b4ed282SShri Abhyankar     PetscReal   wnorm;
1692b4ed282SShri Abhyankar 
1700298fd71SBarry Smith     ierr = VecSetRandom(W,NULL);CHKERRQ(ierr);
1712b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
1722b4ed282SShri Abhyankar     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
1732b4ed282SShri Abhyankar     ierr = MatMult(A,W,work);CHKERRQ(ierr);
1742b4ed282SShri Abhyankar     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
1756bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
1762b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
1774839bfe8SBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
1782b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1792b4ed282SShri Abhyankar   }
1802b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1812b4ed282SShri Abhyankar }
1822b4ed282SShri Abhyankar 
1832b4ed282SShri Abhyankar /*
1842b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
1852b4ed282SShri Abhyankar */
1862b4ed282SShri Abhyankar #undef __FUNCT__
1872b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private"
1882b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
1892b4ed282SShri Abhyankar {
1902b4ed282SShri Abhyankar   PetscReal      a1,a2;
1912b4ed282SShri Abhyankar   PetscErrorCode ierr;
192ace3abfcSBarry Smith   PetscBool      hastranspose;
1932b4ed282SShri Abhyankar 
1942b4ed282SShri Abhyankar   PetscFunctionBegin;
1952b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1962b4ed282SShri Abhyankar   if (hastranspose) {
1972b4ed282SShri Abhyankar     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
1982b4ed282SShri Abhyankar     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
1992b4ed282SShri Abhyankar 
2002b4ed282SShri Abhyankar     /* Compute || J^T W|| */
2012b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
2022b4ed282SShri Abhyankar     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
2032b4ed282SShri Abhyankar     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
2042b4ed282SShri Abhyankar     if (a1 != 0.0) {
2054839bfe8SBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
2062b4ed282SShri Abhyankar     }
2072b4ed282SShri Abhyankar   }
2082b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2092b4ed282SShri Abhyankar }
2102b4ed282SShri Abhyankar 
2112b4ed282SShri Abhyankar /*
2122b4ed282SShri Abhyankar   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
2132b4ed282SShri Abhyankar 
2142b4ed282SShri Abhyankar   Notes:
2152b4ed282SShri Abhyankar   The convergence criterion currently implemented is
2162b4ed282SShri Abhyankar   merit < abstol
2172b4ed282SShri Abhyankar   merit < rtol*merit_initial
2182b4ed282SShri Abhyankar */
2192b4ed282SShri Abhyankar #undef __FUNCT__
2202b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI"
2217fe79bc4SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
2222b4ed282SShri Abhyankar {
2232b4ed282SShri Abhyankar   PetscErrorCode ierr;
2242b4ed282SShri Abhyankar 
2252b4ed282SShri Abhyankar   PetscFunctionBegin;
2262b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2272b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
2282b4ed282SShri Abhyankar 
2292b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
2302b4ed282SShri Abhyankar 
2312b4ed282SShri Abhyankar   if (!it) {
2322b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
2337fe79bc4SShri Abhyankar     snes->ttol = fnorm*snes->rtol;
2342b4ed282SShri Abhyankar   }
2357fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
2362b4ed282SShri Abhyankar     ierr    = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
2372b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
2387fe79bc4SShri Abhyankar   } else if (fnorm < snes->abstol) {
2394839bfe8SBarry Smith     ierr    = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
2402b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
2412b4ed282SShri Abhyankar   } else if (snes->nfuncs >= snes->max_funcs) {
2422b4ed282SShri Abhyankar     ierr    = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
2432b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
2442b4ed282SShri Abhyankar   }
2452b4ed282SShri Abhyankar 
2462b4ed282SShri Abhyankar   if (it && !*reason) {
2477fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2484839bfe8SBarry Smith       ierr    = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
2492b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2502b4ed282SShri Abhyankar     }
2512b4ed282SShri Abhyankar   }
2522b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2532b4ed282SShri Abhyankar }
2542b4ed282SShri Abhyankar 
255ee657d29SShri Abhyankar 
256c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
257c1a5e521SShri Abhyankar /*
258c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
259c1a5e521SShri Abhyankar 
260c1a5e521SShri Abhyankar    Input Parameters:
261c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
262c1a5e521SShri Abhyankar 
263c1a5e521SShri Abhyankar    Output Parameters:
264c1a5e521SShri Abhyankar .  X - Bound projected X
265c1a5e521SShri Abhyankar 
266c1a5e521SShri Abhyankar */
267c1a5e521SShri Abhyankar 
268c1a5e521SShri Abhyankar #undef __FUNCT__
269c1a5e521SShri Abhyankar #define __FUNCT__ "SNESVIProjectOntoBounds"
270c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
271c1a5e521SShri Abhyankar {
272c1a5e521SShri Abhyankar   PetscErrorCode    ierr;
273c1a5e521SShri Abhyankar   const PetscScalar *xl,*xu;
274c1a5e521SShri Abhyankar   PetscScalar       *x;
275c1a5e521SShri Abhyankar   PetscInt          i,n;
276c1a5e521SShri Abhyankar 
277c1a5e521SShri Abhyankar   PetscFunctionBegin;
278c1a5e521SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
279c1a5e521SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
280c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
281c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
282c1a5e521SShri Abhyankar 
283c1a5e521SShri Abhyankar   for (i = 0; i<n; i++) {
28410f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
28510f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
286c1a5e521SShri Abhyankar   }
287c1a5e521SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
288c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
289c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
290c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
291c1a5e521SShri Abhyankar }
292c1a5e521SShri Abhyankar 
29369c03620SShri Abhyankar 
29469c03620SShri Abhyankar #undef __FUNCT__
295693ddf92SShri Abhyankar #define __FUNCT__ "SNESVIGetActiveSetIS"
296693ddf92SShri Abhyankar /*
297693ddf92SShri Abhyankar    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
298693ddf92SShri Abhyankar 
299693ddf92SShri Abhyankar    Input parameter
300693ddf92SShri Abhyankar .  snes - the SNES context
301693ddf92SShri Abhyankar .  X    - the snes solution vector
302693ddf92SShri Abhyankar .  F    - the nonlinear function vector
303693ddf92SShri Abhyankar 
304693ddf92SShri Abhyankar    Output parameter
305693ddf92SShri Abhyankar .  ISact - active set index set
306693ddf92SShri Abhyankar  */
307693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS *ISact)
308d950fb63SShri Abhyankar {
309d950fb63SShri Abhyankar   PetscErrorCode    ierr;
310c2fc9fa9SBarry Smith   Vec               Xl=snes->xl,Xu=snes->xu;
311693ddf92SShri Abhyankar   const PetscScalar *x,*f,*xl,*xu;
312693ddf92SShri Abhyankar   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
313d950fb63SShri Abhyankar 
314d950fb63SShri Abhyankar   PetscFunctionBegin;
315d950fb63SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
316d950fb63SShri Abhyankar   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
317fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
318fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
319fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
320fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
321693ddf92SShri Abhyankar   /* Compute active set size */
322d950fb63SShri Abhyankar   for (i=0; i < nlocal;i++) {
32310f5ae6bSBarry 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++;
324d950fb63SShri Abhyankar   }
325693ddf92SShri Abhyankar 
326d950fb63SShri Abhyankar   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
327d950fb63SShri Abhyankar 
328693ddf92SShri Abhyankar   /* Set active set indices */
329d950fb63SShri Abhyankar   for (i=0; i < nlocal; i++) {
33010f5ae6bSBarry 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;
331d950fb63SShri Abhyankar   }
332d950fb63SShri Abhyankar 
333693ddf92SShri Abhyankar   /* Create active set IS */
334ce94432eSBarry Smith   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes),nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
335d950fb63SShri Abhyankar 
336fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
337fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
338fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
339fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
340d950fb63SShri Abhyankar   PetscFunctionReturn(0);
341d950fb63SShri Abhyankar }
342d950fb63SShri Abhyankar 
343693ddf92SShri Abhyankar #undef __FUNCT__
344693ddf92SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS"
345693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact)
346693ddf92SShri Abhyankar {
347693ddf92SShri Abhyankar   PetscErrorCode ierr;
348077aedafSJed Brown   PetscInt       rstart,rend;
349693ddf92SShri Abhyankar 
350693ddf92SShri Abhyankar   PetscFunctionBegin;
351693ddf92SShri Abhyankar   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
352077aedafSJed Brown   ierr = VecGetOwnershipRange(X,&rstart,&rend);CHKERRQ(ierr);
353077aedafSJed Brown   ierr = ISComplement(*ISact,rstart,rend,ISinact);CHKERRQ(ierr);
354693ddf92SShri Abhyankar   PetscFunctionReturn(0);
355693ddf92SShri Abhyankar }
356693ddf92SShri Abhyankar 
357fe835674SShri Abhyankar #undef __FUNCT__
358fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
35910f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm)
360fe835674SShri Abhyankar {
361fe835674SShri Abhyankar   PetscErrorCode    ierr;
362fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
363fe835674SShri Abhyankar   PetscInt          i,n;
364fe835674SShri Abhyankar   PetscReal         rnorm;
365fe835674SShri Abhyankar 
366fe835674SShri Abhyankar   PetscFunctionBegin;
367fe835674SShri Abhyankar   ierr  = VecGetLocalSize(X,&n);CHKERRQ(ierr);
368c2fc9fa9SBarry Smith   ierr  = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
369c2fc9fa9SBarry Smith   ierr  = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
370fe835674SShri Abhyankar   ierr  = VecGetArrayRead(X,&x);CHKERRQ(ierr);
371fe835674SShri Abhyankar   ierr  = VecGetArrayRead(F,&f);CHKERRQ(ierr);
372fe835674SShri Abhyankar   rnorm = 0.0;
373fe835674SShri Abhyankar   for (i=0; i<n; i++) {
37410f5ae6bSBarry 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]);
3758f918527SKarl Rupp   }
376fe835674SShri Abhyankar   ierr   = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
377c2fc9fa9SBarry Smith   ierr   = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
378c2fc9fa9SBarry Smith   ierr   = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
379fe835674SShri Abhyankar   ierr   = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
380ce94432eSBarry Smith   ierr   = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
38162d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
382fe835674SShri Abhyankar   PetscFunctionReturn(0);
383fe835674SShri Abhyankar }
384fe835674SShri Abhyankar 
38508da532bSDmitry Karpeev #undef __FUNCT__
38608da532bSDmitry Karpeev #define __FUNCT__ "SNESVIDMComputeVariableBounds"
38708da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
38808da532bSDmitry Karpeev {
38908da532bSDmitry Karpeev   PetscErrorCode ierr;
3906e111a19SKarl Rupp 
39108da532bSDmitry Karpeev   PetscFunctionBegin;
39208da532bSDmitry Karpeev   ierr = DMComputeVariableBounds(snes->dm, xl, xu);CHKERRQ(ierr);
39308da532bSDmitry Karpeev   PetscFunctionReturn(0);
39408da532bSDmitry Karpeev }
39508da532bSDmitry Karpeev 
3962f63a38cSShri Abhyankar 
3972b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3982b4ed282SShri Abhyankar /*
399c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
4002b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
4012b4ed282SShri Abhyankar 
4022b4ed282SShri Abhyankar    Input Parameter:
4032b4ed282SShri Abhyankar .  snes - the SNES context
4042b4ed282SShri Abhyankar .  x - the solution vector
4052b4ed282SShri Abhyankar 
4062b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
4072b4ed282SShri Abhyankar 
4082b4ed282SShri Abhyankar    Notes:
4092b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
4102b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
4112b4ed282SShri Abhyankar    the call to SNESSolve().
4122b4ed282SShri Abhyankar  */
4132b4ed282SShri Abhyankar #undef __FUNCT__
4142b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
4152b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
4162b4ed282SShri Abhyankar {
4172b4ed282SShri Abhyankar   PetscErrorCode ierr;
4182b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
4192b4ed282SShri Abhyankar 
4202b4ed282SShri Abhyankar   PetscFunctionBegin;
4219bd66eb0SPeter Brune   ierr = SNESDefaultGetWork(snes,1);CHKERRQ(ierr);
4226cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
4232b4ed282SShri Abhyankar 
42408da532bSDmitry Karpeev   if (!snes->ops->computevariablebounds && snes->dm) {
425a201590fSDmitry Karpeev     PetscBool flag;
426a201590fSDmitry Karpeev     ierr = DMHasVariableBounds(snes->dm, &flag);CHKERRQ(ierr);
4271aa26658SKarl Rupp 
42808da532bSDmitry Karpeev     snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
42908da532bSDmitry Karpeev   }
430a201590fSDmitry Karpeev   if (!snes->usersetbounds) {
431c2fc9fa9SBarry Smith     if (snes->ops->computevariablebounds) {
432c2fc9fa9SBarry Smith       if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
433c2fc9fa9SBarry Smith       if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
434c2fc9fa9SBarry Smith       ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
4351aa26658SKarl Rupp     } else if (!snes->xl && !snes->xu) {
4362176524cSBarry Smith       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
437c2fc9fa9SBarry Smith       ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr);
438c2fc9fa9SBarry Smith       ierr = VecSet(snes->xl,SNES_VI_NINF);CHKERRQ(ierr);
439c2fc9fa9SBarry Smith       ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr);
440c2fc9fa9SBarry Smith       ierr = VecSet(snes->xu,SNES_VI_INF);CHKERRQ(ierr);
4412b4ed282SShri Abhyankar     } else {
4422b4ed282SShri Abhyankar       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
4432b4ed282SShri Abhyankar       ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
444c2fc9fa9SBarry Smith       ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr);
445c2fc9fa9SBarry Smith       ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr);
4462b4ed282SShri 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]))
4472b4ed282SShri 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.");
4482b4ed282SShri Abhyankar     }
449a201590fSDmitry Karpeev   }
4502b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4512b4ed282SShri Abhyankar }
4522b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
4532176524cSBarry Smith #undef __FUNCT__
4542176524cSBarry Smith #define __FUNCT__ "SNESReset_VI"
4552176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes)
4562176524cSBarry Smith {
4572176524cSBarry Smith   PetscErrorCode ierr;
4582176524cSBarry Smith 
4592176524cSBarry Smith   PetscFunctionBegin;
460c2fc9fa9SBarry Smith   ierr                = VecDestroy(&snes->xl);CHKERRQ(ierr);
461c2fc9fa9SBarry Smith   ierr                = VecDestroy(&snes->xu);CHKERRQ(ierr);
4622d6615e8SDmitry Karpeev   snes->usersetbounds = PETSC_FALSE;
4632176524cSBarry Smith   PetscFunctionReturn(0);
4642176524cSBarry Smith }
4652176524cSBarry Smith 
4662b4ed282SShri Abhyankar /*
4672b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
4682b4ed282SShri Abhyankar    with SNESCreate_VI().
4692b4ed282SShri Abhyankar 
4702b4ed282SShri Abhyankar    Input Parameter:
4712b4ed282SShri Abhyankar .  snes - the SNES context
4722b4ed282SShri Abhyankar 
4732b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4742b4ed282SShri Abhyankar  */
4752b4ed282SShri Abhyankar #undef __FUNCT__
4762b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI"
4772b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
4782b4ed282SShri Abhyankar {
4792b4ed282SShri Abhyankar   PetscErrorCode ierr;
4802b4ed282SShri Abhyankar 
4812b4ed282SShri Abhyankar   PetscFunctionBegin;
4822b4ed282SShri Abhyankar   ierr = PetscFree(snes->data);CHKERRQ(ierr);
4832b4ed282SShri Abhyankar 
4842b4ed282SShri Abhyankar   /* clear composed functions */
4850298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",NULL);CHKERRQ(ierr);
4860298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",NULL);CHKERRQ(ierr);
4872b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4882b4ed282SShri Abhyankar }
4897fe79bc4SShri Abhyankar 
4905559b345SBarry Smith #undef __FUNCT__
4915559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds"
4925559b345SBarry Smith /*@
4932b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
4942b4ed282SShri Abhyankar 
4952b4ed282SShri Abhyankar    Input Parameters:
4962b4ed282SShri Abhyankar .  snes - the SNES context.
4972b4ed282SShri Abhyankar .  xl   - lower bound.
4982b4ed282SShri Abhyankar .  xu   - upper bound.
4992b4ed282SShri Abhyankar 
5002b4ed282SShri Abhyankar    Notes:
5012b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
50201f0b76bSBarry Smith    SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp().
50384c105d7SBarry Smith 
5042bd2b0e6SSatish Balay    Level: advanced
5052bd2b0e6SSatish Balay 
5065559b345SBarry Smith @*/
50769c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
5082b4ed282SShri Abhyankar {
50961589011SJed Brown   PetscErrorCode ierr,(*f)(SNES,Vec,Vec);
5102b4ed282SShri Abhyankar 
5112b4ed282SShri Abhyankar   PetscFunctionBegin;
5122b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
5132b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
5142b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
51561589011SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
516f450aa47SBarry Smith   if (!f) {ierr = SNESSetType(snes,SNESVINEWTONRSLS);CHKERRQ(ierr);}
51761589011SJed Brown   ierr                = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr);
518a201590fSDmitry Karpeev   snes->usersetbounds = PETSC_TRUE;
51961589011SJed Brown   PetscFunctionReturn(0);
52061589011SJed Brown }
52161589011SJed Brown 
52261589011SJed Brown EXTERN_C_BEGIN
52361589011SJed Brown #undef __FUNCT__
52461589011SJed Brown #define __FUNCT__ "SNESVISetVariableBounds_VI"
52561589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
52661589011SJed Brown {
52761589011SJed Brown   PetscErrorCode    ierr;
52861589011SJed Brown   const PetscScalar *xxl,*xxu;
52961589011SJed Brown   PetscInt          i,n, cnt = 0;
53061589011SJed Brown 
53161589011SJed Brown   PetscFunctionBegin;
5320298fd71SBarry Smith   ierr = SNESGetFunction(snes,&snes->vec_func,NULL,NULL);CHKERRQ(ierr);
533a63bb30eSJed Brown   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
534077aedafSJed Brown   {
535077aedafSJed Brown     PetscInt xlN,xuN,N;
536077aedafSJed Brown     ierr = VecGetSize(xl,&xlN);CHKERRQ(ierr);
537077aedafSJed Brown     ierr = VecGetSize(xu,&xuN);CHKERRQ(ierr);
538077aedafSJed Brown     ierr = VecGetSize(snes->vec_func,&N);CHKERRQ(ierr);
539077aedafSJed Brown     if (xlN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xlN,N);
540077aedafSJed Brown     if (xuN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xuN,N);
541077aedafSJed Brown   }
542f450aa47SBarry Smith   ierr     = SNESSetType(snes,SNESVINEWTONRSLS);CHKERRQ(ierr);
5432176524cSBarry Smith   ierr     = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
5442176524cSBarry Smith   ierr     = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
545c2fc9fa9SBarry Smith   ierr     = VecDestroy(&snes->xl);CHKERRQ(ierr);
546c2fc9fa9SBarry Smith   ierr     = VecDestroy(&snes->xu);CHKERRQ(ierr);
547c2fc9fa9SBarry Smith   snes->xl = xl;
548c2fc9fa9SBarry Smith   snes->xu = xu;
5496fd67740SBarry Smith   ierr     = VecGetLocalSize(xl,&n);CHKERRQ(ierr);
5506fd67740SBarry Smith   ierr     = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr);
5516fd67740SBarry Smith   ierr     = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr);
5521aa26658SKarl Rupp   for (i=0; i<n; i++) cnt += ((xxl[i] != SNES_VI_NINF) || (xxu[i] != SNES_VI_INF));
5531aa26658SKarl Rupp 
554ce94432eSBarry Smith   ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
5556fd67740SBarry Smith   ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr);
5566fd67740SBarry Smith   ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr);
5572b4ed282SShri Abhyankar   PetscFunctionReturn(0);
5582b4ed282SShri Abhyankar }
55961589011SJed Brown EXTERN_C_END
56092c02d66SPeter Brune 
5612b4ed282SShri Abhyankar #undef __FUNCT__
562c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VI"
563c2fc9fa9SBarry Smith PetscErrorCode SNESSetFromOptions_VI(SNES snes)
5642b4ed282SShri Abhyankar {
5652b4ed282SShri Abhyankar   PetscErrorCode ierr;
566c2fc9fa9SBarry Smith   PetscBool      flg;
567639ea3faSPeter Brune   SNESLineSearch linesearch;
5682b4ed282SShri Abhyankar 
5692b4ed282SShri Abhyankar   PetscFunctionBegin;
570c2fc9fa9SBarry Smith   ierr = PetscOptionsHead("SNES VI options");CHKERRQ(ierr);
571c2fc9fa9SBarry Smith   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
572c2fc9fa9SBarry Smith   if (flg) {
573c2fc9fa9SBarry Smith     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
5742b4ed282SShri Abhyankar   }
575639ea3faSPeter Brune   if (!snes->linesearch) {
576639ea3faSPeter Brune     ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
577639ea3faSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr);
578639ea3faSPeter Brune     ierr = SNESLineSearchBTSetAlpha(linesearch, 0.0);CHKERRQ(ierr);
579639ea3faSPeter Brune   }
580c2fc9fa9SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
581ed70e96aSJungho Lee   PetscFunctionReturn(0);
582ed70e96aSJungho Lee }
583