xref: /petsc/src/snes/impls/vi/vi.c (revision 08da532b2310464c3bf489aef60f0ccd6d61141d)
12b4ed282SShri Abhyankar 
2c2fc9fa9SBarry Smith #include <../include/private/snesimpl.h> /*I "petscsnes.h" I*/
3c6db04a5SJed Brown #include <../include/private/kspimpl.h>
4c6db04a5SJed Brown #include <../include/private/matimpl.h>
5d0af7cd3SBarry Smith #include <../include/private/dmimpl.h>
6d0af7cd3SBarry Smith 
7d0af7cd3SBarry Smith #undef __FUNCT__
82176524cSBarry Smith #define __FUNCT__ "SNESVISetComputeVariableBounds"
92176524cSBarry Smith /*@C
102176524cSBarry Smith    SNESVISetComputeVariableBounds - Sets a function  that is called to compute the variable bounds
112176524cSBarry Smith 
122176524cSBarry Smith    Input parameter
132176524cSBarry Smith +  snes - the SNES context
142176524cSBarry Smith -  compute - computes the bounds
152176524cSBarry Smith 
162bd2b0e6SSatish Balay    Level: advanced
172176524cSBarry Smith 
18aab9d709SJed Brown @*/
1977cdaf69SJed Brown PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
202176524cSBarry Smith {
2161589011SJed Brown   PetscErrorCode   ierr,(*f)(SNES,PetscErrorCode(*)(SNES,Vec,Vec));
222176524cSBarry Smith 
232176524cSBarry Smith   PetscFunctionBegin;
2461589011SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2561589011SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
2661589011SJed Brown   if (!f) {ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);}
2761589011SJed Brown   ierr = PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode(*)(SNES,Vec,Vec)),(snes,compute));CHKERRQ(ierr);
282176524cSBarry Smith   PetscFunctionReturn(0);
292176524cSBarry Smith }
302176524cSBarry Smith 
3161589011SJed Brown EXTERN_C_BEGIN
3261589011SJed Brown #undef __FUNCT__
3361589011SJed Brown #define __FUNCT__ "SNESVISetComputeVariableBounds_VI"
3461589011SJed Brown PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
3561589011SJed Brown {
3661589011SJed Brown   PetscFunctionBegin;
3761589011SJed Brown   snes->ops->computevariablebounds = compute;
3861589011SJed Brown   PetscFunctionReturn(0);
3961589011SJed Brown }
4061589011SJed Brown EXTERN_C_END
412176524cSBarry Smith 
422176524cSBarry Smith #undef __FUNCT__
43ed70e96aSJungho Lee #define __FUNCT__ "SNESVIComputeInactiveSetIS"
44d0af7cd3SBarry Smith /*
45ed70e96aSJungho Lee    SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables
46d0af7cd3SBarry Smith 
47d0af7cd3SBarry Smith    Input parameter
48d0af7cd3SBarry Smith .  snes - the SNES context
49d0af7cd3SBarry Smith .  X    - the snes solution vector
50d0af7cd3SBarry Smith 
51d0af7cd3SBarry Smith    Output parameter
52d0af7cd3SBarry Smith .  ISact - active set index set
53d0af7cd3SBarry Smith 
54d0af7cd3SBarry Smith  */
55ed70e96aSJungho Lee PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS* inact)
56d0af7cd3SBarry Smith {
57d0af7cd3SBarry Smith   PetscErrorCode    ierr;
58d0af7cd3SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
59d0af7cd3SBarry Smith   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
60d0af7cd3SBarry Smith 
61d0af7cd3SBarry Smith   PetscFunctionBegin;
62d0af7cd3SBarry Smith   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
63d0af7cd3SBarry Smith   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
64d0af7cd3SBarry Smith   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
65d0af7cd3SBarry Smith   ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr);
66d0af7cd3SBarry Smith   ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr);
67d0af7cd3SBarry Smith   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
68d0af7cd3SBarry Smith   /* Compute inactive set size */
69d0af7cd3SBarry Smith   for (i=0; i < nlocal;i++) {
70d0af7cd3SBarry 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++;
71d0af7cd3SBarry Smith   }
72d0af7cd3SBarry Smith 
73d0af7cd3SBarry Smith   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
74d0af7cd3SBarry Smith 
75d0af7cd3SBarry Smith   /* Set inactive set indices */
76d0af7cd3SBarry Smith   for(i=0; i < nlocal; i++) {
77d0af7cd3SBarry 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;
78d0af7cd3SBarry Smith   }
79d0af7cd3SBarry Smith 
80d0af7cd3SBarry Smith    /* Create inactive set IS */
81d0af7cd3SBarry Smith   ierr = ISCreateGeneral(((PetscObject)upper)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr);
82d0af7cd3SBarry Smith 
83d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
84d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr);
85d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr);
86d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
87d0af7cd3SBarry Smith   PetscFunctionReturn(0);
88d0af7cd3SBarry Smith }
89d0af7cd3SBarry Smith 
903c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/
912b4ed282SShri Abhyankar 
929308c137SBarry Smith #undef __FUNCT__
939308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI"
947087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
959308c137SBarry Smith {
969308c137SBarry Smith   PetscErrorCode    ierr;
97649052a6SBarry Smith   PetscViewer        viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
989308c137SBarry Smith   const PetscScalar  *x,*xl,*xu,*f;
996fd67740SBarry Smith   PetscInt           i,n,act[2] = {0,0},fact[2],N;
1006a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
1016a9e2d46SJungho Lee   PetscInt           act_bound[2] = {0,0},fact_bound[2];
1029308c137SBarry Smith   PetscReal          rnorm,fnorm;
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]++;
117e400df7aSBarry Smith     else SETERRQ(((PetscObject)snes)->comm,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);
128d9822059SBarry Smith   ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
12921a4c9b1SBarry Smith   ierr = MPI_Allreduce(act,fact,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
1306a9e2d46SJungho Lee   ierr = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
131f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1326fd67740SBarry Smith 
133649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
134c2fc9fa9SBarry Smith   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),((double)(fact[0]+fact[1]))/((double)snes->ntruebounds));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 
1682b4ed282SShri Abhyankar     ierr = VecSetRandom(W,PETSC_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 /*
2102b4ed282SShri Abhyankar   SNESDefaultConverged_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__
2182b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI"
2197fe79bc4SShri Abhyankar PetscErrorCode SNESDefaultConverged_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 
324d950fb63SShri Abhyankar   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&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 */
33262298a1eSBarry Smith   ierr = ISCreateGeneral(((PetscObject)snes)->comm,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;
346693ddf92SShri Abhyankar 
347693ddf92SShri Abhyankar   PetscFunctionBegin;
348693ddf92SShri Abhyankar   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
349693ddf92SShri Abhyankar   ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr);
350693ddf92SShri Abhyankar   PetscFunctionReturn(0);
351693ddf92SShri Abhyankar }
352693ddf92SShri Abhyankar 
353fe835674SShri Abhyankar #undef __FUNCT__
354fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
35510f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X,PetscReal *fnorm)
356fe835674SShri Abhyankar {
357fe835674SShri Abhyankar   PetscErrorCode    ierr;
358fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
359fe835674SShri Abhyankar   PetscInt          i,n;
360fe835674SShri Abhyankar   PetscReal         rnorm;
361fe835674SShri Abhyankar 
362fe835674SShri Abhyankar   PetscFunctionBegin;
363fe835674SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
364c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
365c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
366fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
367fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
368fe835674SShri Abhyankar   rnorm = 0.0;
369fe835674SShri Abhyankar   for (i=0; i<n; i++) {
37010f5ae6bSBarry 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]);
371fe835674SShri Abhyankar   }
372fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
373c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
374c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
375fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
376d9822059SBarry Smith   ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
37762d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
378fe835674SShri Abhyankar   PetscFunctionReturn(0);
379fe835674SShri Abhyankar }
380fe835674SShri Abhyankar 
381*08da532bSDmitry Karpeev #undef __FUNCT__
382*08da532bSDmitry Karpeev #define __FUNCT__ "SNESVIDMComputeVariableBounds"
383*08da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
384*08da532bSDmitry Karpeev {
385*08da532bSDmitry Karpeev   PetscErrorCode     ierr;
386*08da532bSDmitry Karpeev   PetscFunctionBegin;
387*08da532bSDmitry Karpeev   ierr = DMComputeVariableBounds(snes->dm, xl, xu); CHKERRQ(ierr);
388*08da532bSDmitry Karpeev   PetscFunctionReturn(0);
389*08da532bSDmitry Karpeev }
390*08da532bSDmitry Karpeev 
3912f63a38cSShri Abhyankar 
3922b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3932b4ed282SShri Abhyankar /*
394c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3952b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3962b4ed282SShri Abhyankar 
3972b4ed282SShri Abhyankar    Input Parameter:
3982b4ed282SShri Abhyankar .  snes - the SNES context
3992b4ed282SShri Abhyankar .  x - the solution vector
4002b4ed282SShri Abhyankar 
4012b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
4022b4ed282SShri Abhyankar 
4032b4ed282SShri Abhyankar    Notes:
4042b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
4052b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
4062b4ed282SShri Abhyankar    the call to SNESSolve().
4072b4ed282SShri Abhyankar  */
4082b4ed282SShri Abhyankar #undef __FUNCT__
4092b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
4102b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
4112b4ed282SShri Abhyankar {
4122b4ed282SShri Abhyankar   PetscErrorCode ierr;
4132b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
4142b4ed282SShri Abhyankar 
4152b4ed282SShri Abhyankar   PetscFunctionBegin;
41658c9b817SLisandro Dalcin 
41758c9b817SLisandro Dalcin   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
4182b4ed282SShri Abhyankar 
419*08da532bSDmitry Karpeev   if(!snes->ops->computevariablebounds && snes->dm) {
420*08da532bSDmitry Karpeev     snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
421*08da532bSDmitry 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 
4802b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
481c2fc9fa9SBarry Smith 
482c2fc9fa9SBarry Smith /*
483c2fc9fa9SBarry Smith 
484c2fc9fa9SBarry Smith       These line searches are common for all the VI solvers
485c2fc9fa9SBarry Smith */
486c2fc9fa9SBarry Smith extern PetscErrorCode SNESSolve_VISS(SNES);
487c2fc9fa9SBarry Smith 
4882b4ed282SShri Abhyankar #undef __FUNCT__
4892b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI"
4902b4ed282SShri Abhyankar 
4912b4ed282SShri Abhyankar /*
4927fe79bc4SShri Abhyankar   This routine does not actually do a line search but it takes a full newton
4937fe79bc4SShri Abhyankar   step while ensuring that the new iterates remain within the constraints.
4942b4ed282SShri Abhyankar 
4952b4ed282SShri Abhyankar */
4961e633543SBarry Smith PetscErrorCode SNESLineSearchNo_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
4972b4ed282SShri Abhyankar {
4982b4ed282SShri Abhyankar   PetscErrorCode ierr;
499ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5002b4ed282SShri Abhyankar 
5012b4ed282SShri Abhyankar   PetscFunctionBegin;
5022b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
5032b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5042b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
5052b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
506c1a5e521SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
5070661309fSPeter Brune   if (snes->ops->postcheckstep) {
5080661309fSPeter Brune    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
509c1a5e521SShri Abhyankar   }
510c1a5e521SShri Abhyankar   if (changed_y) {
511c1a5e521SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
5127fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
513c1a5e521SShri Abhyankar   }
514c1a5e521SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
515c1a5e521SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
516c1a5e521SShri Abhyankar   if (!snes->domainerror) {
517c2fc9fa9SBarry Smith     if (snes->ops->solve != SNESSolve_VISS) {
5187fe79bc4SShri Abhyankar        ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
5197fe79bc4SShri Abhyankar     } else {
520c1a5e521SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
5217fe79bc4SShri Abhyankar     }
52262cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
523c1a5e521SShri Abhyankar   }
5240661309fSPeter Brune   if (snes->ls_monitor) {
5250661309fSPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
5260661309fSPeter Brune     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
5270661309fSPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
528c1a5e521SShri Abhyankar   }
529c1a5e521SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
530c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
531c1a5e521SShri Abhyankar }
532c1a5e521SShri Abhyankar 
533c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
534c1a5e521SShri Abhyankar #undef __FUNCT__
5352b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI"
5362b4ed282SShri Abhyankar 
5372b4ed282SShri Abhyankar /*
5382b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
5392b4ed282SShri Abhyankar */
5401e633543SBarry Smith PetscErrorCode SNESLineSearchNoNorms_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
5412b4ed282SShri Abhyankar {
5422b4ed282SShri Abhyankar   PetscErrorCode ierr;
543ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5442b4ed282SShri Abhyankar 
5452b4ed282SShri Abhyankar   PetscFunctionBegin;
5462b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
5472b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5482b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
5497fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
5500661309fSPeter Brune   if (snes->ops->postcheckstep) {
5510661309fSPeter Brune    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
5522b4ed282SShri Abhyankar   }
5532b4ed282SShri Abhyankar   if (changed_y) {
5542b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
5557fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
5562b4ed282SShri Abhyankar   }
5572b4ed282SShri Abhyankar 
5582b4ed282SShri Abhyankar   /* don't evaluate function the last time through */
5592b4ed282SShri Abhyankar   if (snes->iter < snes->max_its-1) {
5602b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
5612b4ed282SShri Abhyankar   }
5622b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5632b4ed282SShri Abhyankar   PetscFunctionReturn(0);
5642b4ed282SShri Abhyankar }
5657fe79bc4SShri Abhyankar 
5662b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
567c2fc9fa9SBarry Smith 
5682b4ed282SShri Abhyankar #undef __FUNCT__
5692b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI"
5702b4ed282SShri Abhyankar /*
5717fe79bc4SShri Abhyankar   This routine implements a cubic line search while doing a projection on the variable bounds
5722b4ed282SShri Abhyankar */
5731e633543SBarry Smith PetscErrorCode SNESLineSearchCubic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
5742b4ed282SShri Abhyankar {
575fe835674SShri Abhyankar   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
576fe835674SShri Abhyankar   PetscReal      minlambda,lambda,lambdatemp;
577fe835674SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
578fe835674SShri Abhyankar   PetscScalar    cinitslope;
579fe835674SShri Abhyankar #endif
580fe835674SShri Abhyankar   PetscErrorCode ierr;
581fe835674SShri Abhyankar   PetscInt       count;
582fe835674SShri Abhyankar   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
583fe835674SShri Abhyankar   MPI_Comm       comm;
584fe835674SShri Abhyankar 
585fe835674SShri Abhyankar   PetscFunctionBegin;
586fe835674SShri Abhyankar   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
587fe835674SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
588fe835674SShri Abhyankar   *flag   = PETSC_TRUE;
589fe835674SShri Abhyankar 
590fe835674SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
591fe835674SShri Abhyankar   if (*ynorm == 0.0) {
5920661309fSPeter Brune     if (snes->ls_monitor) {
5930661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
5940661309fSPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
5950661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
596fe835674SShri Abhyankar     }
597fe835674SShri Abhyankar     *gnorm = fnorm;
598fe835674SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
599fe835674SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
600fe835674SShri Abhyankar     *flag  = PETSC_FALSE;
601fe835674SShri Abhyankar     goto theend1;
602fe835674SShri Abhyankar   }
6030661309fSPeter Brune   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
6040661309fSPeter Brune     if (snes->ls_monitor) {
6050661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
60662d1f40fSBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Scaling step by %g old ynorm %g\n",(double)(snes->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr);
6070661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
608fe835674SShri Abhyankar     }
6090661309fSPeter Brune     ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
6100661309fSPeter Brune     *ynorm = snes->maxstep;
611fe835674SShri Abhyankar   }
612fe835674SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
6130661309fSPeter Brune   minlambda = snes->steptol/rellength;
614fe835674SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
615fe835674SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
616fe835674SShri Abhyankar   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
617fe835674SShri Abhyankar   initslope = PetscRealPart(cinitslope);
618fe835674SShri Abhyankar #else
619fe835674SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
620fe835674SShri Abhyankar #endif
621fe835674SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
622fe835674SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
623fe835674SShri Abhyankar 
624fe835674SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
625fe835674SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
626fe835674SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
627fe835674SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
628fe835674SShri Abhyankar     *flag = PETSC_FALSE;
629fe835674SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
630fe835674SShri Abhyankar     goto theend1;
631fe835674SShri Abhyankar   }
632fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
633c2fc9fa9SBarry Smith   if (snes->ops->solve != SNESSolve_VISS) {
634fe835674SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
6357fe79bc4SShri Abhyankar   } else {
6367fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
6377fe79bc4SShri Abhyankar   }
638fe835674SShri Abhyankar   if (snes->domainerror) {
639fe835674SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
640fe835674SShri Abhyankar     PetscFunctionReturn(0);
641fe835674SShri Abhyankar   }
64262cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6430661309fSPeter Brune   ierr = PetscInfo4(snes,"Initial fnorm %g gnorm %g alpha %g initslope %g\n",(double)fnorm,(double)*gnorm,(double)snes->ls_alpha,(double)initslope);CHKERRQ(ierr);
6440661309fSPeter Brune   if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* Sufficient reduction */
6450661309fSPeter Brune     if (snes->ls_monitor) {
6460661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
64762d1f40fSBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr);
6480661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
649fe835674SShri Abhyankar     }
650fe835674SShri Abhyankar     goto theend1;
651fe835674SShri Abhyankar   }
652fe835674SShri Abhyankar 
653fe835674SShri Abhyankar   /* Fit points with quadratic */
654fe835674SShri Abhyankar   lambda     = 1.0;
655fe835674SShri Abhyankar   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
656fe835674SShri Abhyankar   lambdaprev = lambda;
657fe835674SShri Abhyankar   gnormprev  = *gnorm;
658fe835674SShri Abhyankar   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
659fe835674SShri Abhyankar   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
660fe835674SShri Abhyankar   else                         lambda = lambdatemp;
661fe835674SShri Abhyankar 
662fe835674SShri Abhyankar   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
663fe835674SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
664fe835674SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
665fe835674SShri Abhyankar     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
666fe835674SShri Abhyankar     *flag = PETSC_FALSE;
667fe835674SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
668fe835674SShri Abhyankar     goto theend1;
669fe835674SShri Abhyankar   }
670fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
671c2fc9fa9SBarry Smith   if (snes->ops->solve != SNESSolve_VISS) {
672fe835674SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
6737fe79bc4SShri Abhyankar   } else {
6747fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
6757fe79bc4SShri Abhyankar   }
676fe835674SShri Abhyankar   if (snes->domainerror) {
677fe835674SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
678fe835674SShri Abhyankar     PetscFunctionReturn(0);
679fe835674SShri Abhyankar   }
68062cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6810661309fSPeter Brune   if (snes->ls_monitor) {
6820661309fSPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
68362d1f40fSBarry Smith     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: gnorm after quadratic fit %g\n",(double)(*gnorm));CHKERRQ(ierr);
6840661309fSPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
685fe835674SShri Abhyankar   }
6860661309fSPeter Brune   if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* sufficient reduction */
6870661309fSPeter Brune     if (snes->ls_monitor) {
6880661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
6894839bfe8SBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
6900661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
691fe835674SShri Abhyankar     }
692fe835674SShri Abhyankar     goto theend1;
693fe835674SShri Abhyankar   }
694fe835674SShri Abhyankar 
695fe835674SShri Abhyankar   /* Fit points with cubic */
696fe835674SShri Abhyankar   count = 1;
697fe835674SShri Abhyankar   while (PETSC_TRUE) {
698fe835674SShri Abhyankar     if (lambda <= minlambda) {
6990661309fSPeter Brune       if (snes->ls_monitor) {
7000661309fSPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
7010661309fSPeter Brune  	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
70262d1f40fSBarry Smith 	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)minlambda,(double)lambda,(double)initslope);CHKERRQ(ierr);
7030661309fSPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
704fe835674SShri Abhyankar       }
705fe835674SShri Abhyankar       *flag = PETSC_FALSE;
706fe835674SShri Abhyankar       break;
707fe835674SShri Abhyankar     }
708fe835674SShri Abhyankar     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
709fe835674SShri Abhyankar     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
710fe835674SShri Abhyankar     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
711fe835674SShri Abhyankar     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
712fe835674SShri Abhyankar     d  = b*b - 3*a*initslope;
713fe835674SShri Abhyankar     if (d < 0.0) d = 0.0;
714fe835674SShri Abhyankar     if (a == 0.0) {
715fe835674SShri Abhyankar       lambdatemp = -initslope/(2.0*b);
716fe835674SShri Abhyankar     } else {
71762d1f40fSBarry Smith       lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
718fe835674SShri Abhyankar     }
719fe835674SShri Abhyankar     lambdaprev = lambda;
720fe835674SShri Abhyankar     gnormprev  = *gnorm;
721fe835674SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
722fe835674SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
723fe835674SShri Abhyankar     else                         lambda     = lambdatemp;
724fe835674SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
725fe835674SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
726fe835674SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
727fe835674SShri Abhyankar       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
7284839bfe8SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)lambda,(double)initslope);CHKERRQ(ierr);
729fe835674SShri Abhyankar       *flag = PETSC_FALSE;
730fe835674SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
731fe835674SShri Abhyankar       break;
732fe835674SShri Abhyankar     }
733fe835674SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
734c2fc9fa9SBarry Smith     if (snes->ops->solve != SNESSolve_VISS) {
735fe835674SShri Abhyankar       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
7367fe79bc4SShri Abhyankar     } else {
7377fe79bc4SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7387fe79bc4SShri Abhyankar     }
739fe835674SShri Abhyankar     if (snes->domainerror) {
740fe835674SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
741fe835674SShri Abhyankar       PetscFunctionReturn(0);
742fe835674SShri Abhyankar     }
74362cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7440661309fSPeter Brune     if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* is reduction enough? */
7450661309fSPeter Brune       if (snes->ls_monitor) {
74662d1f40fSBarry Smith 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %g lambda=%18.16e\n",(double)(*gnorm),(double)lambda);CHKERRQ(ierr);
747fe835674SShri Abhyankar       }
748fe835674SShri Abhyankar       break;
749fe835674SShri Abhyankar     } else {
7500661309fSPeter Brune       if (snes->ls_monitor) {
75162d1f40fSBarry Smith         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %g lambda=%18.16e\n",(double)(*gnorm),(double)lambda);CHKERRQ(ierr);
752fe835674SShri Abhyankar       }
753fe835674SShri Abhyankar     }
754fe835674SShri Abhyankar     count++;
755fe835674SShri Abhyankar   }
756fe835674SShri Abhyankar   theend1:
757fe835674SShri Abhyankar   /* Optional user-defined check for line search step validity */
7580661309fSPeter Brune   if (snes->ops->postcheckstep && *flag) {
7590661309fSPeter Brune     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
760fe835674SShri Abhyankar     if (changed_y) {
761fe835674SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
762fe835674SShri Abhyankar       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
763fe835674SShri Abhyankar     }
764fe835674SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
765fe835674SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
766c2fc9fa9SBarry Smith       if (snes->ops->solve != SNESSolve_VISS) {
767fe835674SShri Abhyankar         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
7687fe79bc4SShri Abhyankar       } else {
7697fe79bc4SShri Abhyankar         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7707fe79bc4SShri Abhyankar       }
771fe835674SShri Abhyankar       if (snes->domainerror) {
772fe835674SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
773fe835674SShri Abhyankar         PetscFunctionReturn(0);
774fe835674SShri Abhyankar       }
77562cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
776fe835674SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
777fe835674SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
778fe835674SShri Abhyankar     }
779fe835674SShri Abhyankar   }
780fe835674SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
781fe835674SShri Abhyankar   PetscFunctionReturn(0);
782fe835674SShri Abhyankar }
783fe835674SShri Abhyankar 
7842b4ed282SShri Abhyankar #undef __FUNCT__
7852b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI"
7862b4ed282SShri Abhyankar /*
7877fe79bc4SShri Abhyankar   This routine does a quadratic line search while keeping the iterates within the variable bounds
7882b4ed282SShri Abhyankar */
7891e633543SBarry Smith PetscErrorCode SNESLineSearchQuadratic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag)
7902b4ed282SShri Abhyankar {
7912b4ed282SShri Abhyankar   /*
7922b4ed282SShri Abhyankar      Note that for line search purposes we work with with the related
7932b4ed282SShri Abhyankar      minimization problem:
7942b4ed282SShri Abhyankar         min  z(x):  R^n -> R,
7952b4ed282SShri Abhyankar      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
7962b4ed282SShri Abhyankar    */
7972b4ed282SShri Abhyankar   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
7982b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
7992b4ed282SShri Abhyankar   PetscScalar    cinitslope;
8002b4ed282SShri Abhyankar #endif
8012b4ed282SShri Abhyankar   PetscErrorCode ierr;
8022b4ed282SShri Abhyankar   PetscInt       count;
803ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
8042b4ed282SShri Abhyankar 
8052b4ed282SShri Abhyankar   PetscFunctionBegin;
8062b4ed282SShri Abhyankar   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8072b4ed282SShri Abhyankar   *flag   = PETSC_TRUE;
8082b4ed282SShri Abhyankar 
8092b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
8102b4ed282SShri Abhyankar   if (*ynorm == 0.0) {
8110661309fSPeter Brune     if (snes->ls_monitor) {
8120661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8130661309fSPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
8140661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8152b4ed282SShri Abhyankar     }
8162b4ed282SShri Abhyankar     *gnorm = fnorm;
8172b4ed282SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
8182b4ed282SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
8192b4ed282SShri Abhyankar     *flag  = PETSC_FALSE;
8202b4ed282SShri Abhyankar     goto theend2;
8212b4ed282SShri Abhyankar   }
8220661309fSPeter Brune   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
8230661309fSPeter Brune     ierr   = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
8240661309fSPeter Brune     *ynorm = snes->maxstep;
8252b4ed282SShri Abhyankar   }
8262b4ed282SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
8270661309fSPeter Brune   minlambda = snes->steptol/rellength;
8287fe79bc4SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
8292b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
8307fe79bc4SShri Abhyankar   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
8312b4ed282SShri Abhyankar   initslope = PetscRealPart(cinitslope);
8322b4ed282SShri Abhyankar #else
8337fe79bc4SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
8342b4ed282SShri Abhyankar #endif
8352b4ed282SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
8362b4ed282SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
8374839bfe8SBarry Smith   ierr = PetscInfo1(snes,"Initslope %g \n",(double)initslope);CHKERRQ(ierr);
8382b4ed282SShri Abhyankar 
8392b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8407fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
8412b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
8422b4ed282SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
8432b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
8442b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
8452b4ed282SShri Abhyankar     goto theend2;
8462b4ed282SShri Abhyankar   }
8472b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
848c2fc9fa9SBarry Smith   if (snes->ops->solve != SNESSolve_VISS) {
8497fe79bc4SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
8507fe79bc4SShri Abhyankar   } else {
8517fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
8527fe79bc4SShri Abhyankar   }
8532b4ed282SShri Abhyankar   if (snes->domainerror) {
8542b4ed282SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8552b4ed282SShri Abhyankar     PetscFunctionReturn(0);
8562b4ed282SShri Abhyankar   }
85762cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8580661309fSPeter Brune   if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* Sufficient reduction */
8590661309fSPeter Brune     if (snes->ls_monitor) {
8600661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
86162d1f40fSBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr);
8620661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8632b4ed282SShri Abhyankar     }
8642b4ed282SShri Abhyankar     goto theend2;
8652b4ed282SShri Abhyankar   }
8662b4ed282SShri Abhyankar 
8672b4ed282SShri Abhyankar   /* Fit points with quadratic */
8682b4ed282SShri Abhyankar   lambda = 1.0;
8692b4ed282SShri Abhyankar   count = 1;
8702b4ed282SShri Abhyankar   while (PETSC_TRUE) {
8712b4ed282SShri Abhyankar     if (lambda <= minlambda) { /* bad luck; use full step */
8720661309fSPeter Brune       if (snes->ls_monitor) {
8730661309fSPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8740661309fSPeter Brune         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
87562d1f40fSBarry Smith         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)lambda,(double)initslope);CHKERRQ(ierr);
8760661309fSPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8772b4ed282SShri Abhyankar       }
8782b4ed282SShri Abhyankar       ierr = VecCopy(x,w);CHKERRQ(ierr);
8792b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
8802b4ed282SShri Abhyankar       break;
8812b4ed282SShri Abhyankar     }
8822b4ed282SShri Abhyankar     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
8832b4ed282SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
8842b4ed282SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
8852b4ed282SShri Abhyankar     else                         lambda     = lambdatemp;
8862b4ed282SShri Abhyankar 
8872b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
8887fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
8892b4ed282SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
8902b4ed282SShri Abhyankar       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
8914839bfe8SBarry Smith       ierr  = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)lambda,(double)initslope);CHKERRQ(ierr);
8922b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
8932b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
8942b4ed282SShri Abhyankar       break;
8952b4ed282SShri Abhyankar     }
8962b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
8972b4ed282SShri Abhyankar     if (snes->domainerror) {
8982b4ed282SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8992b4ed282SShri Abhyankar       PetscFunctionReturn(0);
9002b4ed282SShri Abhyankar     }
901c2fc9fa9SBarry Smith     if (snes->ops->solve != SNESSolve_VISS) {
9027fe79bc4SShri Abhyankar       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
9037fe79bc4SShri Abhyankar     } else {
9042b4ed282SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
9057fe79bc4SShri Abhyankar     }
90662cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
9070661309fSPeter Brune     if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* sufficient reduction */
9080661309fSPeter Brune       if (snes->ls_monitor) {
9090661309fSPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
91062d1f40fSBarry Smith         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Quadratically determined step, lambda=%g\n",(double)lambda);CHKERRQ(ierr);
9110661309fSPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
9122b4ed282SShri Abhyankar       }
9132b4ed282SShri Abhyankar       break;
9142b4ed282SShri Abhyankar     }
9152b4ed282SShri Abhyankar     count++;
9162b4ed282SShri Abhyankar   }
9172b4ed282SShri Abhyankar   theend2:
9182b4ed282SShri Abhyankar   /* Optional user-defined check for line search step validity */
9190661309fSPeter Brune   if (snes->ops->postcheckstep) {
9200661309fSPeter Brune     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
9212b4ed282SShri Abhyankar     if (changed_y) {
9222b4ed282SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
9237fe79bc4SShri Abhyankar       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
9242b4ed282SShri Abhyankar     }
9252b4ed282SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
9262b4ed282SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);
9272b4ed282SShri Abhyankar       if (snes->domainerror) {
9282b4ed282SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
9292b4ed282SShri Abhyankar         PetscFunctionReturn(0);
9302b4ed282SShri Abhyankar       }
931c2fc9fa9SBarry Smith       if (snes->ops->solve != SNESSolve_VISS) {
9327fe79bc4SShri Abhyankar         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
9337fe79bc4SShri Abhyankar       } else {
9347fe79bc4SShri Abhyankar         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
9357fe79bc4SShri Abhyankar       }
9367fe79bc4SShri Abhyankar 
9372b4ed282SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
9382b4ed282SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
93962cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
9402b4ed282SShri Abhyankar     }
9412b4ed282SShri Abhyankar   }
9422b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
9432b4ed282SShri Abhyankar   PetscFunctionReturn(0);
9442b4ed282SShri Abhyankar }
9452b4ed282SShri Abhyankar 
9462b4ed282SShri Abhyankar 
9475559b345SBarry Smith #undef __FUNCT__
9485559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds"
9495559b345SBarry Smith /*@
9502b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
9512b4ed282SShri Abhyankar 
9522b4ed282SShri Abhyankar    Input Parameters:
9532b4ed282SShri Abhyankar .  snes - the SNES context.
9542b4ed282SShri Abhyankar .  xl   - lower bound.
9552b4ed282SShri Abhyankar .  xu   - upper bound.
9562b4ed282SShri Abhyankar 
9572b4ed282SShri Abhyankar    Notes:
9582b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
95901f0b76bSBarry Smith    SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp().
96084c105d7SBarry Smith 
9612bd2b0e6SSatish Balay    Level: advanced
9622bd2b0e6SSatish Balay 
9635559b345SBarry Smith @*/
96469c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
9652b4ed282SShri Abhyankar {
96661589011SJed Brown   PetscErrorCode   ierr,(*f)(SNES,Vec,Vec);
9672b4ed282SShri Abhyankar 
9682b4ed282SShri Abhyankar   PetscFunctionBegin;
9692b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
9702b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
9712b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
97261589011SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
97361589011SJed Brown   if (!f) {ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);}
97461589011SJed Brown   ierr = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr);
97561589011SJed Brown   PetscFunctionReturn(0);
97661589011SJed Brown }
97761589011SJed Brown 
97861589011SJed Brown EXTERN_C_BEGIN
97961589011SJed Brown #undef __FUNCT__
98061589011SJed Brown #define __FUNCT__ "SNESVISetVariableBounds_VI"
98161589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
98261589011SJed Brown {
98361589011SJed Brown   PetscErrorCode    ierr;
98461589011SJed Brown   const PetscScalar *xxl,*xxu;
98561589011SJed Brown   PetscInt          i,n, cnt = 0;
98661589011SJed Brown 
98761589011SJed Brown   PetscFunctionBegin;
988a63bb30eSJed Brown   ierr = SNESGetFunction(snes,&snes->vec_func,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
989a63bb30eSJed Brown   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
9902b4ed282SShri 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);
9912b4ed282SShri 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);
992c2fc9fa9SBarry Smith   ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);
9932176524cSBarry Smith   ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
9942176524cSBarry Smith   ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
995c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
996c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
997c2fc9fa9SBarry Smith   snes->xl = xl;
998c2fc9fa9SBarry Smith   snes->xu = xu;
9996fd67740SBarry Smith   ierr = VecGetLocalSize(xl,&n);CHKERRQ(ierr);
10006fd67740SBarry Smith   ierr = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr);
10016fd67740SBarry Smith   ierr = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr);
10026fd67740SBarry Smith   for (i=0; i<n; i++) {
10036fd67740SBarry Smith     cnt += ((xxl[i] != SNES_VI_NINF) || (xxu[i] != SNES_VI_INF));
10046fd67740SBarry Smith   }
1005c2fc9fa9SBarry Smith   ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
10066fd67740SBarry Smith   ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr);
10076fd67740SBarry Smith   ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr);
10082b4ed282SShri Abhyankar   PetscFunctionReturn(0);
10092b4ed282SShri Abhyankar }
101061589011SJed Brown EXTERN_C_END
101192c02d66SPeter Brune 
101292c02d66SPeter Brune EXTERN_C_BEGIN
101392c02d66SPeter Brune #undef __FUNCT__
101492c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_VI"
101592c02d66SPeter Brune PetscErrorCode  SNESLineSearchSetType_VI(SNES snes, SNESLineSearchType type)
101692c02d66SPeter Brune {
101792c02d66SPeter Brune   PetscErrorCode ierr;
101892c02d66SPeter Brune   PetscFunctionBegin;
101992c02d66SPeter Brune 
102092c02d66SPeter Brune   switch (type) {
102192c02d66SPeter Brune   case SNES_LS_BASIC:
102292c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
102392c02d66SPeter Brune     break;
102492c02d66SPeter Brune   case SNES_LS_BASIC_NONORMS:
102592c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
102692c02d66SPeter Brune     break;
102792c02d66SPeter Brune   case SNES_LS_QUADRATIC:
102892c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
102992c02d66SPeter Brune     break;
103092c02d66SPeter Brune   case SNES_LS_CUBIC:
103192c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
103292c02d66SPeter Brune     break;
103392c02d66SPeter Brune   default:
103492c02d66SPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
103592c02d66SPeter Brune     break;
103692c02d66SPeter Brune   }
103792c02d66SPeter Brune   snes->ls_type = type;
103892c02d66SPeter Brune   PetscFunctionReturn(0);
103992c02d66SPeter Brune }
104092c02d66SPeter Brune EXTERN_C_END
104192c02d66SPeter Brune 
10422b4ed282SShri Abhyankar #undef __FUNCT__
1043c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VI"
1044c2fc9fa9SBarry Smith PetscErrorCode SNESSetFromOptions_VI(SNES snes)
10452b4ed282SShri Abhyankar {
10462b4ed282SShri Abhyankar   PetscErrorCode  ierr;
1047c2fc9fa9SBarry Smith   PetscBool       flg;
10482b4ed282SShri Abhyankar 
10492b4ed282SShri Abhyankar   PetscFunctionBegin;
1050c2fc9fa9SBarry Smith   ierr = PetscOptionsHead("SNES VI options");CHKERRQ(ierr);
1051c2fc9fa9SBarry Smith   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
1052c2fc9fa9SBarry Smith   if (flg) {
1053c2fc9fa9SBarry Smith     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
10542b4ed282SShri Abhyankar   }
1055c2fc9fa9SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1056ed70e96aSJungho Lee   PetscFunctionReturn(0);
1057ed70e96aSJungho Lee }
1058