xref: /petsc/src/snes/impls/vi/vi.c (revision 9bd66eb0dc31efbcccaa04c7c4d5de114a77b297)
1*9bd66eb0SPeter Brune #include "viimpl.h"  /*I "petscsnes.h" I*/
2d0af7cd3SBarry Smith 
3d0af7cd3SBarry Smith #undef __FUNCT__
42176524cSBarry Smith #define __FUNCT__ "SNESVISetComputeVariableBounds"
52176524cSBarry Smith /*@C
62176524cSBarry Smith    SNESVISetComputeVariableBounds - Sets a function  that is called to compute the variable bounds
72176524cSBarry Smith 
82176524cSBarry Smith    Input parameter
92176524cSBarry Smith +  snes - the SNES context
102176524cSBarry Smith -  compute - computes the bounds
112176524cSBarry Smith 
122bd2b0e6SSatish Balay    Level: advanced
132176524cSBarry Smith 
14aab9d709SJed Brown @*/
1577cdaf69SJed Brown PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
162176524cSBarry Smith {
1761589011SJed Brown   PetscErrorCode   ierr,(*f)(SNES,PetscErrorCode(*)(SNES,Vec,Vec));
182176524cSBarry Smith 
192176524cSBarry Smith   PetscFunctionBegin;
2061589011SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2161589011SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
2261589011SJed Brown   if (!f) {ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);}
2361589011SJed Brown   ierr = PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode(*)(SNES,Vec,Vec)),(snes,compute));CHKERRQ(ierr);
242176524cSBarry Smith   PetscFunctionReturn(0);
252176524cSBarry Smith }
262176524cSBarry Smith 
2761589011SJed Brown EXTERN_C_BEGIN
2861589011SJed Brown #undef __FUNCT__
2961589011SJed Brown #define __FUNCT__ "SNESVISetComputeVariableBounds_VI"
3061589011SJed Brown PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
3161589011SJed Brown {
3261589011SJed Brown   PetscFunctionBegin;
3361589011SJed Brown   snes->ops->computevariablebounds = compute;
3461589011SJed Brown   PetscFunctionReturn(0);
3561589011SJed Brown }
3661589011SJed Brown EXTERN_C_END
372176524cSBarry Smith 
382176524cSBarry Smith #undef __FUNCT__
39ed70e96aSJungho Lee #define __FUNCT__ "SNESVIComputeInactiveSetIS"
40d0af7cd3SBarry Smith /*
41ed70e96aSJungho Lee    SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables
42d0af7cd3SBarry Smith 
43d0af7cd3SBarry Smith    Input parameter
44d0af7cd3SBarry Smith .  snes - the SNES context
45d0af7cd3SBarry Smith .  X    - the snes solution vector
46d0af7cd3SBarry Smith 
47d0af7cd3SBarry Smith    Output parameter
48d0af7cd3SBarry Smith .  ISact - active set index set
49d0af7cd3SBarry Smith 
50d0af7cd3SBarry Smith  */
51ed70e96aSJungho Lee PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS* inact)
52d0af7cd3SBarry Smith {
53d0af7cd3SBarry Smith   PetscErrorCode    ierr;
54d0af7cd3SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
55d0af7cd3SBarry Smith   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
56d0af7cd3SBarry Smith 
57d0af7cd3SBarry Smith   PetscFunctionBegin;
58d0af7cd3SBarry Smith   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
59d0af7cd3SBarry Smith   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
60d0af7cd3SBarry Smith   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
61d0af7cd3SBarry Smith   ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr);
62d0af7cd3SBarry Smith   ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr);
63d0af7cd3SBarry Smith   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
64d0af7cd3SBarry Smith   /* Compute inactive set size */
65d0af7cd3SBarry Smith   for (i=0; i < nlocal;i++) {
66d0af7cd3SBarry Smith     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) nloc_isact++;
67d0af7cd3SBarry Smith   }
68d0af7cd3SBarry Smith 
69d0af7cd3SBarry Smith   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
70d0af7cd3SBarry Smith 
71d0af7cd3SBarry Smith   /* Set inactive set indices */
72d0af7cd3SBarry Smith   for(i=0; i < nlocal; i++) {
73d0af7cd3SBarry Smith     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) idx_act[i1++] = ilow+i;
74d0af7cd3SBarry Smith   }
75d0af7cd3SBarry Smith 
76d0af7cd3SBarry Smith    /* Create inactive set IS */
77d0af7cd3SBarry Smith   ierr = ISCreateGeneral(((PetscObject)upper)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr);
78d0af7cd3SBarry Smith 
79d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
80d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr);
81d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr);
82d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
83d0af7cd3SBarry Smith   PetscFunctionReturn(0);
84d0af7cd3SBarry Smith }
85d0af7cd3SBarry Smith 
863c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/
872b4ed282SShri Abhyankar 
889308c137SBarry Smith #undef __FUNCT__
899308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI"
907087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
919308c137SBarry Smith {
929308c137SBarry Smith   PetscErrorCode    ierr;
93649052a6SBarry Smith   PetscViewer        viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
949308c137SBarry Smith   const PetscScalar  *x,*xl,*xu,*f;
956fd67740SBarry Smith   PetscInt           i,n,act[2] = {0,0},fact[2],N;
966a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
976a9e2d46SJungho Lee   PetscInt           act_bound[2] = {0,0},fact_bound[2];
989308c137SBarry Smith   PetscReal          rnorm,fnorm;
999308c137SBarry Smith 
1009308c137SBarry Smith   PetscFunctionBegin;
1019308c137SBarry Smith   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
1026fd67740SBarry Smith   ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr);
103c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
104c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
1059308c137SBarry Smith   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
1063f731af5SBarry Smith   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
1079308c137SBarry Smith 
1089308c137SBarry Smith   rnorm = 0.0;
1099308c137SBarry Smith   for (i=0; i<n; i++) {
11010f5ae6bSBarry 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]);
111e400df7aSBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++;
112e400df7aSBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++;
113e400df7aSBarry Smith     else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Can never get here");
1149308c137SBarry Smith   }
1156a9e2d46SJungho Lee 
1166a9e2d46SJungho Lee   for (i=0; i<n; i++) {
1176a9e2d46SJungho Lee     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++;
1186a9e2d46SJungho Lee     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++;
1196a9e2d46SJungho Lee   }
1203f731af5SBarry Smith   ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
121c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
122c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
1239308c137SBarry Smith   ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
124d9822059SBarry Smith   ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
12521a4c9b1SBarry Smith   ierr = MPI_Allreduce(act,fact,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
1266a9e2d46SJungho Lee   ierr = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
127f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1286fd67740SBarry Smith 
129649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
130c2fc9fa9SBarry 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);
1316a9e2d46SJungho Lee 
132649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1339308c137SBarry Smith   PetscFunctionReturn(0);
1349308c137SBarry Smith }
1359308c137SBarry Smith 
1362b4ed282SShri Abhyankar /*
1372b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1382b4ed282SShri 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
1392b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1402b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1412b4ed282SShri Abhyankar */
1422b4ed282SShri Abhyankar #undef __FUNCT__
1432b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private"
144ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
1452b4ed282SShri Abhyankar {
1462b4ed282SShri Abhyankar   PetscReal      a1;
1472b4ed282SShri Abhyankar   PetscErrorCode ierr;
148ace3abfcSBarry Smith   PetscBool     hastranspose;
1492b4ed282SShri Abhyankar 
1502b4ed282SShri Abhyankar   PetscFunctionBegin;
1512b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1522b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1532b4ed282SShri Abhyankar   if (hastranspose) {
1542b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1552b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
1562b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
1574839bfe8SBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
1582b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
1592b4ed282SShri Abhyankar   } else {
1602b4ed282SShri Abhyankar     Vec         work;
1612b4ed282SShri Abhyankar     PetscScalar result;
1622b4ed282SShri Abhyankar     PetscReal   wnorm;
1632b4ed282SShri Abhyankar 
1642b4ed282SShri Abhyankar     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
1652b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
1662b4ed282SShri Abhyankar     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
1672b4ed282SShri Abhyankar     ierr = MatMult(A,W,work);CHKERRQ(ierr);
1682b4ed282SShri Abhyankar     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
1696bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
1702b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
1714839bfe8SBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
1722b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1732b4ed282SShri Abhyankar   }
1742b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1752b4ed282SShri Abhyankar }
1762b4ed282SShri Abhyankar 
1772b4ed282SShri Abhyankar /*
1782b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
1792b4ed282SShri Abhyankar */
1802b4ed282SShri Abhyankar #undef __FUNCT__
1812b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private"
1822b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
1832b4ed282SShri Abhyankar {
1842b4ed282SShri Abhyankar   PetscReal      a1,a2;
1852b4ed282SShri Abhyankar   PetscErrorCode ierr;
186ace3abfcSBarry Smith   PetscBool     hastranspose;
1872b4ed282SShri Abhyankar 
1882b4ed282SShri Abhyankar   PetscFunctionBegin;
1892b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1902b4ed282SShri Abhyankar   if (hastranspose) {
1912b4ed282SShri Abhyankar     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
1922b4ed282SShri Abhyankar     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
1932b4ed282SShri Abhyankar 
1942b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1952b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
1962b4ed282SShri Abhyankar     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
1972b4ed282SShri Abhyankar     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
1982b4ed282SShri Abhyankar     if (a1 != 0.0) {
1994839bfe8SBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
2002b4ed282SShri Abhyankar     }
2012b4ed282SShri Abhyankar   }
2022b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2032b4ed282SShri Abhyankar }
2042b4ed282SShri Abhyankar 
2052b4ed282SShri Abhyankar /*
2062b4ed282SShri Abhyankar   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
2072b4ed282SShri Abhyankar 
2082b4ed282SShri Abhyankar   Notes:
2092b4ed282SShri Abhyankar   The convergence criterion currently implemented is
2102b4ed282SShri Abhyankar   merit < abstol
2112b4ed282SShri Abhyankar   merit < rtol*merit_initial
2122b4ed282SShri Abhyankar */
2132b4ed282SShri Abhyankar #undef __FUNCT__
2142b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI"
2157fe79bc4SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
2162b4ed282SShri Abhyankar {
2172b4ed282SShri Abhyankar   PetscErrorCode ierr;
2182b4ed282SShri Abhyankar 
2192b4ed282SShri Abhyankar   PetscFunctionBegin;
2202b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2212b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
2222b4ed282SShri Abhyankar 
2232b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
2242b4ed282SShri Abhyankar 
2252b4ed282SShri Abhyankar   if (!it) {
2262b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
2277fe79bc4SShri Abhyankar     snes->ttol = fnorm*snes->rtol;
2282b4ed282SShri Abhyankar   }
2297fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
2302b4ed282SShri Abhyankar     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
2312b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
2327fe79bc4SShri Abhyankar   } else if (fnorm < snes->abstol) {
2334839bfe8SBarry Smith     ierr = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
2342b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
2352b4ed282SShri Abhyankar   } else if (snes->nfuncs >= snes->max_funcs) {
2362b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
2372b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
2382b4ed282SShri Abhyankar   }
2392b4ed282SShri Abhyankar 
2402b4ed282SShri Abhyankar   if (it && !*reason) {
2417fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2424839bfe8SBarry Smith       ierr = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
2432b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2442b4ed282SShri Abhyankar     }
2452b4ed282SShri Abhyankar   }
2462b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2472b4ed282SShri Abhyankar }
2482b4ed282SShri Abhyankar 
249ee657d29SShri Abhyankar 
250c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
251c1a5e521SShri Abhyankar /*
252c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
253c1a5e521SShri Abhyankar 
254c1a5e521SShri Abhyankar    Input Parameters:
255c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
256c1a5e521SShri Abhyankar 
257c1a5e521SShri Abhyankar    Output Parameters:
258c1a5e521SShri Abhyankar .  X - Bound projected X
259c1a5e521SShri Abhyankar 
260c1a5e521SShri Abhyankar */
261c1a5e521SShri Abhyankar 
262c1a5e521SShri Abhyankar #undef __FUNCT__
263c1a5e521SShri Abhyankar #define __FUNCT__ "SNESVIProjectOntoBounds"
264c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
265c1a5e521SShri Abhyankar {
266c1a5e521SShri Abhyankar   PetscErrorCode    ierr;
267c1a5e521SShri Abhyankar   const PetscScalar *xl,*xu;
268c1a5e521SShri Abhyankar   PetscScalar       *x;
269c1a5e521SShri Abhyankar   PetscInt          i,n;
270c1a5e521SShri Abhyankar 
271c1a5e521SShri Abhyankar   PetscFunctionBegin;
272c1a5e521SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
273c1a5e521SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
274c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
275c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
276c1a5e521SShri Abhyankar 
277c1a5e521SShri Abhyankar   for(i = 0;i<n;i++) {
27810f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
27910f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
280c1a5e521SShri Abhyankar   }
281c1a5e521SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
282c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
283c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
284c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
285c1a5e521SShri Abhyankar }
286c1a5e521SShri Abhyankar 
28769c03620SShri Abhyankar 
28869c03620SShri Abhyankar #undef __FUNCT__
289693ddf92SShri Abhyankar #define __FUNCT__ "SNESVIGetActiveSetIS"
290693ddf92SShri Abhyankar /*
291693ddf92SShri Abhyankar    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
292693ddf92SShri Abhyankar 
293693ddf92SShri Abhyankar    Input parameter
294693ddf92SShri Abhyankar .  snes - the SNES context
295693ddf92SShri Abhyankar .  X    - the snes solution vector
296693ddf92SShri Abhyankar .  F    - the nonlinear function vector
297693ddf92SShri Abhyankar 
298693ddf92SShri Abhyankar    Output parameter
299693ddf92SShri Abhyankar .  ISact - active set index set
300693ddf92SShri Abhyankar  */
301693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS* ISact)
302d950fb63SShri Abhyankar {
303d950fb63SShri Abhyankar   PetscErrorCode   ierr;
304c2fc9fa9SBarry Smith   Vec               Xl=snes->xl,Xu=snes->xu;
305693ddf92SShri Abhyankar   const PetscScalar *x,*f,*xl,*xu;
306693ddf92SShri Abhyankar   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
307d950fb63SShri Abhyankar 
308d950fb63SShri Abhyankar   PetscFunctionBegin;
309d950fb63SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
310d950fb63SShri Abhyankar   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
311fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
312fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
313fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
314fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
315693ddf92SShri Abhyankar   /* Compute active set size */
316d950fb63SShri Abhyankar   for (i=0; i < nlocal;i++) {
31710f5ae6bSBarry 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++;
318d950fb63SShri Abhyankar   }
319693ddf92SShri Abhyankar 
320d950fb63SShri Abhyankar   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
321d950fb63SShri Abhyankar 
322693ddf92SShri Abhyankar   /* Set active set indices */
323d950fb63SShri Abhyankar   for(i=0; i < nlocal; i++) {
32410f5ae6bSBarry 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;
325d950fb63SShri Abhyankar   }
326d950fb63SShri Abhyankar 
327693ddf92SShri Abhyankar    /* Create active set IS */
32862298a1eSBarry Smith   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
329d950fb63SShri Abhyankar 
330fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
331fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
332fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
333fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
334d950fb63SShri Abhyankar   PetscFunctionReturn(0);
335d950fb63SShri Abhyankar }
336d950fb63SShri Abhyankar 
337693ddf92SShri Abhyankar #undef __FUNCT__
338693ddf92SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS"
339693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS* ISact,IS* ISinact)
340693ddf92SShri Abhyankar {
341693ddf92SShri Abhyankar   PetscErrorCode     ierr;
342693ddf92SShri Abhyankar 
343693ddf92SShri Abhyankar   PetscFunctionBegin;
344693ddf92SShri Abhyankar   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
345693ddf92SShri Abhyankar   ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr);
346693ddf92SShri Abhyankar   PetscFunctionReturn(0);
347693ddf92SShri Abhyankar }
348693ddf92SShri Abhyankar 
349fe835674SShri Abhyankar #undef __FUNCT__
350fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
35110f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm)
352fe835674SShri Abhyankar {
353fe835674SShri Abhyankar   PetscErrorCode    ierr;
354fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
355fe835674SShri Abhyankar   PetscInt          i,n;
356fe835674SShri Abhyankar   PetscReal         rnorm;
357fe835674SShri Abhyankar 
358fe835674SShri Abhyankar   PetscFunctionBegin;
359fe835674SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
360c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
361c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
362fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
363fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
364fe835674SShri Abhyankar   rnorm = 0.0;
365fe835674SShri Abhyankar   for (i=0; i<n; i++) {
36610f5ae6bSBarry 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]);
367fe835674SShri Abhyankar   }
368fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
369c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
370c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
371fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
372d9822059SBarry Smith   ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
37362d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
374fe835674SShri Abhyankar   PetscFunctionReturn(0);
375fe835674SShri Abhyankar }
376fe835674SShri Abhyankar 
37708da532bSDmitry Karpeev #undef __FUNCT__
37808da532bSDmitry Karpeev #define __FUNCT__ "SNESVIDMComputeVariableBounds"
37908da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
38008da532bSDmitry Karpeev {
38108da532bSDmitry Karpeev   PetscErrorCode     ierr;
38208da532bSDmitry Karpeev   PetscFunctionBegin;
38308da532bSDmitry Karpeev   ierr = DMComputeVariableBounds(snes->dm, xl, xu); CHKERRQ(ierr);
38408da532bSDmitry Karpeev   PetscFunctionReturn(0);
38508da532bSDmitry Karpeev }
38608da532bSDmitry Karpeev 
3872f63a38cSShri Abhyankar 
3882b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3892b4ed282SShri Abhyankar /*
390c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3912b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3922b4ed282SShri Abhyankar 
3932b4ed282SShri Abhyankar    Input Parameter:
3942b4ed282SShri Abhyankar .  snes - the SNES context
3952b4ed282SShri Abhyankar .  x - the solution vector
3962b4ed282SShri Abhyankar 
3972b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
3982b4ed282SShri Abhyankar 
3992b4ed282SShri Abhyankar    Notes:
4002b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
4012b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
4022b4ed282SShri Abhyankar    the call to SNESSolve().
4032b4ed282SShri Abhyankar  */
4042b4ed282SShri Abhyankar #undef __FUNCT__
4052b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
4062b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
4072b4ed282SShri Abhyankar {
4082b4ed282SShri Abhyankar   PetscErrorCode ierr;
4092b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
4102b4ed282SShri Abhyankar 
4112b4ed282SShri Abhyankar   PetscFunctionBegin;
41258c9b817SLisandro Dalcin 
413*9bd66eb0SPeter Brune   ierr = SNESDefaultGetWork(snes,1);CHKERRQ(ierr);
4146cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
4152b4ed282SShri Abhyankar 
41608da532bSDmitry Karpeev   if(!snes->ops->computevariablebounds && snes->dm) {
41708da532bSDmitry Karpeev     snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
41808da532bSDmitry Karpeev   }
419c2fc9fa9SBarry Smith   if (snes->ops->computevariablebounds) {
420c2fc9fa9SBarry Smith     if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
421c2fc9fa9SBarry Smith     if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
422c2fc9fa9SBarry Smith     ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
423c2fc9fa9SBarry Smith   } else if (!snes->xl && !snes->xu) {
4242176524cSBarry Smith     /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
425c2fc9fa9SBarry Smith     ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr);
426c2fc9fa9SBarry Smith     ierr = VecSet(snes->xl,SNES_VI_NINF);CHKERRQ(ierr);
427c2fc9fa9SBarry Smith     ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr);
428c2fc9fa9SBarry Smith     ierr = VecSet(snes->xu,SNES_VI_INF);CHKERRQ(ierr);
4292b4ed282SShri Abhyankar   } else {
4302b4ed282SShri Abhyankar     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
4312b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
432c2fc9fa9SBarry Smith     ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr);
433c2fc9fa9SBarry Smith     ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr);
4342b4ed282SShri 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]))
4352b4ed282SShri 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.");
4362b4ed282SShri Abhyankar   }
437abcba341SShri Abhyankar 
4382b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4392b4ed282SShri Abhyankar }
4402b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
4412176524cSBarry Smith #undef __FUNCT__
4422176524cSBarry Smith #define __FUNCT__ "SNESReset_VI"
4432176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes)
4442176524cSBarry Smith {
4452176524cSBarry Smith   PetscErrorCode ierr;
4462176524cSBarry Smith 
4472176524cSBarry Smith   PetscFunctionBegin;
448c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
449c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
4502176524cSBarry Smith   PetscFunctionReturn(0);
4512176524cSBarry Smith }
4522176524cSBarry Smith 
4532b4ed282SShri Abhyankar /*
4542b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
4552b4ed282SShri Abhyankar    with SNESCreate_VI().
4562b4ed282SShri Abhyankar 
4572b4ed282SShri Abhyankar    Input Parameter:
4582b4ed282SShri Abhyankar .  snes - the SNES context
4592b4ed282SShri Abhyankar 
4602b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4612b4ed282SShri Abhyankar  */
4622b4ed282SShri Abhyankar #undef __FUNCT__
4632b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI"
4642b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
4652b4ed282SShri Abhyankar {
4662b4ed282SShri Abhyankar   PetscErrorCode ierr;
4672b4ed282SShri Abhyankar 
4682b4ed282SShri Abhyankar   PetscFunctionBegin;
4692b4ed282SShri Abhyankar   ierr = PetscFree(snes->data);CHKERRQ(ierr);
4702b4ed282SShri Abhyankar 
4712b4ed282SShri Abhyankar   /* clear composed functions */
4722b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
473c92abb8aSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
4742b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4752b4ed282SShri Abhyankar }
4767fe79bc4SShri Abhyankar 
4772b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
478c2fc9fa9SBarry Smith 
479c2fc9fa9SBarry Smith /*
480c2fc9fa9SBarry Smith 
481c2fc9fa9SBarry Smith       These line searches are common for all the VI solvers
482c2fc9fa9SBarry Smith */
483c2fc9fa9SBarry Smith extern PetscErrorCode SNESSolve_VISS(SNES);
484c2fc9fa9SBarry Smith 
4852b4ed282SShri Abhyankar #undef __FUNCT__
4862b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI"
4872b4ed282SShri Abhyankar 
4882b4ed282SShri Abhyankar /*
4897fe79bc4SShri Abhyankar   This routine does not actually do a line search but it takes a full newton
4907fe79bc4SShri Abhyankar   step while ensuring that the new iterates remain within the constraints.
4912b4ed282SShri Abhyankar 
4922b4ed282SShri Abhyankar */
4931e633543SBarry 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)
4942b4ed282SShri Abhyankar {
4952b4ed282SShri Abhyankar   PetscErrorCode ierr;
496ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
4972b4ed282SShri Abhyankar 
4982b4ed282SShri Abhyankar   PetscFunctionBegin;
4992b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
5002b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5012b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
5022b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
503c1a5e521SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
5040661309fSPeter Brune   if (snes->ops->postcheckstep) {
5050661309fSPeter Brune    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
506c1a5e521SShri Abhyankar   }
507c1a5e521SShri Abhyankar   if (changed_y) {
508c1a5e521SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
5097fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
510c1a5e521SShri Abhyankar   }
511c1a5e521SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
512c1a5e521SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
513c1a5e521SShri Abhyankar   if (!snes->domainerror) {
514c2fc9fa9SBarry Smith     if (snes->ops->solve != SNESSolve_VISS) {
5157fe79bc4SShri Abhyankar        ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
5167fe79bc4SShri Abhyankar     } else {
517c1a5e521SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
5187fe79bc4SShri Abhyankar     }
51962cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
520c1a5e521SShri Abhyankar   }
5210661309fSPeter Brune   if (snes->ls_monitor) {
5220661309fSPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
5230661309fSPeter Brune     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
5240661309fSPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
525c1a5e521SShri Abhyankar   }
526c1a5e521SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
527c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
528c1a5e521SShri Abhyankar }
529c1a5e521SShri Abhyankar 
530c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
531c1a5e521SShri Abhyankar #undef __FUNCT__
5322b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI"
5332b4ed282SShri Abhyankar 
5342b4ed282SShri Abhyankar /*
5352b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
5362b4ed282SShri Abhyankar */
5371e633543SBarry 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)
5382b4ed282SShri Abhyankar {
5392b4ed282SShri Abhyankar   PetscErrorCode ierr;
540ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5412b4ed282SShri Abhyankar 
5422b4ed282SShri Abhyankar   PetscFunctionBegin;
5432b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
5442b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5452b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
5467fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
5470661309fSPeter Brune   if (snes->ops->postcheckstep) {
5480661309fSPeter Brune    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
5492b4ed282SShri Abhyankar   }
5502b4ed282SShri Abhyankar   if (changed_y) {
5512b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
5527fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
5532b4ed282SShri Abhyankar   }
5542b4ed282SShri Abhyankar 
5552b4ed282SShri Abhyankar   /* don't evaluate function the last time through */
5562b4ed282SShri Abhyankar   if (snes->iter < snes->max_its-1) {
5572b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
5582b4ed282SShri Abhyankar   }
5592b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5602b4ed282SShri Abhyankar   PetscFunctionReturn(0);
5612b4ed282SShri Abhyankar }
5627fe79bc4SShri Abhyankar 
5632b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
564c2fc9fa9SBarry Smith 
5652b4ed282SShri Abhyankar #undef __FUNCT__
5662b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI"
5672b4ed282SShri Abhyankar /*
5687fe79bc4SShri Abhyankar   This routine implements a cubic line search while doing a projection on the variable bounds
5692b4ed282SShri Abhyankar */
5701e633543SBarry 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)
5712b4ed282SShri Abhyankar {
572fe835674SShri Abhyankar   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
573fe835674SShri Abhyankar   PetscReal      minlambda,lambda,lambdatemp;
574fe835674SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
575fe835674SShri Abhyankar   PetscScalar    cinitslope;
576fe835674SShri Abhyankar #endif
577fe835674SShri Abhyankar   PetscErrorCode ierr;
578fe835674SShri Abhyankar   PetscInt       count;
579fe835674SShri Abhyankar   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
580fe835674SShri Abhyankar   MPI_Comm       comm;
581fe835674SShri Abhyankar 
582fe835674SShri Abhyankar   PetscFunctionBegin;
583fe835674SShri Abhyankar   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
584fe835674SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
585fe835674SShri Abhyankar   *flag   = PETSC_TRUE;
586fe835674SShri Abhyankar 
587fe835674SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
588fe835674SShri Abhyankar   if (*ynorm == 0.0) {
5890661309fSPeter Brune     if (snes->ls_monitor) {
5900661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
5910661309fSPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
5920661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
593fe835674SShri Abhyankar     }
594fe835674SShri Abhyankar     *gnorm = fnorm;
595fe835674SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
596fe835674SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
597fe835674SShri Abhyankar     *flag  = PETSC_FALSE;
598fe835674SShri Abhyankar     goto theend1;
599fe835674SShri Abhyankar   }
6000661309fSPeter Brune   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
6010661309fSPeter Brune     if (snes->ls_monitor) {
6020661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
60362d1f40fSBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Scaling step by %g old ynorm %g\n",(double)(snes->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr);
6040661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
605fe835674SShri Abhyankar     }
6060661309fSPeter Brune     ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
6070661309fSPeter Brune     *ynorm = snes->maxstep;
608fe835674SShri Abhyankar   }
609fe835674SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
6100661309fSPeter Brune   minlambda = snes->steptol/rellength;
611fe835674SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
612fe835674SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
613fe835674SShri Abhyankar   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
614fe835674SShri Abhyankar   initslope = PetscRealPart(cinitslope);
615fe835674SShri Abhyankar #else
616fe835674SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
617fe835674SShri Abhyankar #endif
618fe835674SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
619fe835674SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
620fe835674SShri Abhyankar 
621fe835674SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
622fe835674SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
623fe835674SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
624fe835674SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
625fe835674SShri Abhyankar     *flag = PETSC_FALSE;
626fe835674SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
627fe835674SShri Abhyankar     goto theend1;
628fe835674SShri Abhyankar   }
629fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
630c2fc9fa9SBarry Smith   if (snes->ops->solve != SNESSolve_VISS) {
631fe835674SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
6327fe79bc4SShri Abhyankar   } else {
6337fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
6347fe79bc4SShri Abhyankar   }
635fe835674SShri Abhyankar   if (snes->domainerror) {
636fe835674SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
637fe835674SShri Abhyankar     PetscFunctionReturn(0);
638fe835674SShri Abhyankar   }
63962cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6400661309fSPeter 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);
6410661309fSPeter Brune   if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* Sufficient reduction */
6420661309fSPeter Brune     if (snes->ls_monitor) {
6430661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
64462d1f40fSBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr);
6450661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
646fe835674SShri Abhyankar     }
647fe835674SShri Abhyankar     goto theend1;
648fe835674SShri Abhyankar   }
649fe835674SShri Abhyankar 
650fe835674SShri Abhyankar   /* Fit points with quadratic */
651fe835674SShri Abhyankar   lambda     = 1.0;
652fe835674SShri Abhyankar   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
653fe835674SShri Abhyankar   lambdaprev = lambda;
654fe835674SShri Abhyankar   gnormprev  = *gnorm;
655fe835674SShri Abhyankar   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
656fe835674SShri Abhyankar   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
657fe835674SShri Abhyankar   else                         lambda = lambdatemp;
658fe835674SShri Abhyankar 
659fe835674SShri Abhyankar   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
660fe835674SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
661fe835674SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
662fe835674SShri Abhyankar     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
663fe835674SShri Abhyankar     *flag = PETSC_FALSE;
664fe835674SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
665fe835674SShri Abhyankar     goto theend1;
666fe835674SShri Abhyankar   }
667fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
668c2fc9fa9SBarry Smith   if (snes->ops->solve != SNESSolve_VISS) {
669fe835674SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
6707fe79bc4SShri Abhyankar   } else {
6717fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
6727fe79bc4SShri Abhyankar   }
673fe835674SShri Abhyankar   if (snes->domainerror) {
674fe835674SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
675fe835674SShri Abhyankar     PetscFunctionReturn(0);
676fe835674SShri Abhyankar   }
67762cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6780661309fSPeter Brune   if (snes->ls_monitor) {
6790661309fSPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
68062d1f40fSBarry Smith     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: gnorm after quadratic fit %g\n",(double)(*gnorm));CHKERRQ(ierr);
6810661309fSPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
682fe835674SShri Abhyankar   }
6830661309fSPeter Brune   if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* sufficient reduction */
6840661309fSPeter Brune     if (snes->ls_monitor) {
6850661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
6864839bfe8SBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
6870661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
688fe835674SShri Abhyankar     }
689fe835674SShri Abhyankar     goto theend1;
690fe835674SShri Abhyankar   }
691fe835674SShri Abhyankar 
692fe835674SShri Abhyankar   /* Fit points with cubic */
693fe835674SShri Abhyankar   count = 1;
694fe835674SShri Abhyankar   while (PETSC_TRUE) {
695fe835674SShri Abhyankar     if (lambda <= minlambda) {
6960661309fSPeter Brune       if (snes->ls_monitor) {
6970661309fSPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
6980661309fSPeter Brune  	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
69962d1f40fSBarry 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);
7000661309fSPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
701fe835674SShri Abhyankar       }
702fe835674SShri Abhyankar       *flag = PETSC_FALSE;
703fe835674SShri Abhyankar       break;
704fe835674SShri Abhyankar     }
705fe835674SShri Abhyankar     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
706fe835674SShri Abhyankar     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
707fe835674SShri Abhyankar     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
708fe835674SShri Abhyankar     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
709fe835674SShri Abhyankar     d  = b*b - 3*a*initslope;
710fe835674SShri Abhyankar     if (d < 0.0) d = 0.0;
711fe835674SShri Abhyankar     if (a == 0.0) {
712fe835674SShri Abhyankar       lambdatemp = -initslope/(2.0*b);
713fe835674SShri Abhyankar     } else {
71462d1f40fSBarry Smith       lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
715fe835674SShri Abhyankar     }
716fe835674SShri Abhyankar     lambdaprev = lambda;
717fe835674SShri Abhyankar     gnormprev  = *gnorm;
718fe835674SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
719fe835674SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
720fe835674SShri Abhyankar     else                         lambda     = lambdatemp;
721fe835674SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
722fe835674SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
723fe835674SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
724fe835674SShri Abhyankar       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
7254839bfe8SBarry 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);
726fe835674SShri Abhyankar       *flag = PETSC_FALSE;
727fe835674SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
728fe835674SShri Abhyankar       break;
729fe835674SShri Abhyankar     }
730fe835674SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
731c2fc9fa9SBarry Smith     if (snes->ops->solve != SNESSolve_VISS) {
732fe835674SShri Abhyankar       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
7337fe79bc4SShri Abhyankar     } else {
7347fe79bc4SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7357fe79bc4SShri Abhyankar     }
736fe835674SShri Abhyankar     if (snes->domainerror) {
737fe835674SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
738fe835674SShri Abhyankar       PetscFunctionReturn(0);
739fe835674SShri Abhyankar     }
74062cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7410661309fSPeter Brune     if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* is reduction enough? */
7420661309fSPeter Brune       if (snes->ls_monitor) {
74362d1f40fSBarry Smith 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %g lambda=%18.16e\n",(double)(*gnorm),(double)lambda);CHKERRQ(ierr);
744fe835674SShri Abhyankar       }
745fe835674SShri Abhyankar       break;
746fe835674SShri Abhyankar     } else {
7470661309fSPeter Brune       if (snes->ls_monitor) {
74862d1f40fSBarry 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);
749fe835674SShri Abhyankar       }
750fe835674SShri Abhyankar     }
751fe835674SShri Abhyankar     count++;
752fe835674SShri Abhyankar   }
753fe835674SShri Abhyankar   theend1:
754fe835674SShri Abhyankar   /* Optional user-defined check for line search step validity */
7550661309fSPeter Brune   if (snes->ops->postcheckstep && *flag) {
7560661309fSPeter Brune     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
757fe835674SShri Abhyankar     if (changed_y) {
758fe835674SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
759fe835674SShri Abhyankar       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
760fe835674SShri Abhyankar     }
761fe835674SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
762fe835674SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
763c2fc9fa9SBarry Smith       if (snes->ops->solve != SNESSolve_VISS) {
764fe835674SShri Abhyankar         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
7657fe79bc4SShri Abhyankar       } else {
7667fe79bc4SShri Abhyankar         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7677fe79bc4SShri Abhyankar       }
768fe835674SShri Abhyankar       if (snes->domainerror) {
769fe835674SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
770fe835674SShri Abhyankar         PetscFunctionReturn(0);
771fe835674SShri Abhyankar       }
77262cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
773fe835674SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
774fe835674SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
775fe835674SShri Abhyankar     }
776fe835674SShri Abhyankar   }
777fe835674SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
778fe835674SShri Abhyankar   PetscFunctionReturn(0);
779fe835674SShri Abhyankar }
780fe835674SShri Abhyankar 
7812b4ed282SShri Abhyankar #undef __FUNCT__
7822b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI"
7832b4ed282SShri Abhyankar /*
7847fe79bc4SShri Abhyankar   This routine does a quadratic line search while keeping the iterates within the variable bounds
7852b4ed282SShri Abhyankar */
7861e633543SBarry 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)
7872b4ed282SShri Abhyankar {
7882b4ed282SShri Abhyankar   /*
7892b4ed282SShri Abhyankar      Note that for line search purposes we work with with the related
7902b4ed282SShri Abhyankar      minimization problem:
7912b4ed282SShri Abhyankar         min  z(x):  R^n -> R,
7922b4ed282SShri Abhyankar      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
7932b4ed282SShri Abhyankar    */
7942b4ed282SShri Abhyankar   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
7952b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
7962b4ed282SShri Abhyankar   PetscScalar    cinitslope;
7972b4ed282SShri Abhyankar #endif
7982b4ed282SShri Abhyankar   PetscErrorCode ierr;
7992b4ed282SShri Abhyankar   PetscInt       count;
800ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
8012b4ed282SShri Abhyankar 
8022b4ed282SShri Abhyankar   PetscFunctionBegin;
8032b4ed282SShri Abhyankar   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8042b4ed282SShri Abhyankar   *flag   = PETSC_TRUE;
8052b4ed282SShri Abhyankar 
8062b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
8072b4ed282SShri Abhyankar   if (*ynorm == 0.0) {
8080661309fSPeter Brune     if (snes->ls_monitor) {
8090661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8100661309fSPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
8110661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8122b4ed282SShri Abhyankar     }
8132b4ed282SShri Abhyankar     *gnorm = fnorm;
8142b4ed282SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
8152b4ed282SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
8162b4ed282SShri Abhyankar     *flag  = PETSC_FALSE;
8172b4ed282SShri Abhyankar     goto theend2;
8182b4ed282SShri Abhyankar   }
8190661309fSPeter Brune   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
8200661309fSPeter Brune     ierr   = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
8210661309fSPeter Brune     *ynorm = snes->maxstep;
8222b4ed282SShri Abhyankar   }
8232b4ed282SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
8240661309fSPeter Brune   minlambda = snes->steptol/rellength;
8257fe79bc4SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
8262b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
8277fe79bc4SShri Abhyankar   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
8282b4ed282SShri Abhyankar   initslope = PetscRealPart(cinitslope);
8292b4ed282SShri Abhyankar #else
8307fe79bc4SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
8312b4ed282SShri Abhyankar #endif
8322b4ed282SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
8332b4ed282SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
8344839bfe8SBarry Smith   ierr = PetscInfo1(snes,"Initslope %g \n",(double)initslope);CHKERRQ(ierr);
8352b4ed282SShri Abhyankar 
8362b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8377fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
8382b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
8392b4ed282SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
8402b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
8412b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
8422b4ed282SShri Abhyankar     goto theend2;
8432b4ed282SShri Abhyankar   }
8442b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
845c2fc9fa9SBarry Smith   if (snes->ops->solve != SNESSolve_VISS) {
8467fe79bc4SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
8477fe79bc4SShri Abhyankar   } else {
8487fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
8497fe79bc4SShri Abhyankar   }
8502b4ed282SShri Abhyankar   if (snes->domainerror) {
8512b4ed282SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8522b4ed282SShri Abhyankar     PetscFunctionReturn(0);
8532b4ed282SShri Abhyankar   }
85462cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8550661309fSPeter Brune   if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* Sufficient reduction */
8560661309fSPeter Brune     if (snes->ls_monitor) {
8570661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
85862d1f40fSBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr);
8590661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8602b4ed282SShri Abhyankar     }
8612b4ed282SShri Abhyankar     goto theend2;
8622b4ed282SShri Abhyankar   }
8632b4ed282SShri Abhyankar 
8642b4ed282SShri Abhyankar   /* Fit points with quadratic */
8652b4ed282SShri Abhyankar   lambda = 1.0;
8662b4ed282SShri Abhyankar   count = 1;
8672b4ed282SShri Abhyankar   while (PETSC_TRUE) {
8682b4ed282SShri Abhyankar     if (lambda <= minlambda) { /* bad luck; use full step */
8690661309fSPeter Brune       if (snes->ls_monitor) {
8700661309fSPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8710661309fSPeter Brune         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
87262d1f40fSBarry 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);
8730661309fSPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8742b4ed282SShri Abhyankar       }
8752b4ed282SShri Abhyankar       ierr = VecCopy(x,w);CHKERRQ(ierr);
8762b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
8772b4ed282SShri Abhyankar       break;
8782b4ed282SShri Abhyankar     }
8792b4ed282SShri Abhyankar     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
8802b4ed282SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
8812b4ed282SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
8822b4ed282SShri Abhyankar     else                         lambda     = lambdatemp;
8832b4ed282SShri Abhyankar 
8842b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
8857fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
8862b4ed282SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
8872b4ed282SShri Abhyankar       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
8884839bfe8SBarry 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);
8892b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
8902b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
8912b4ed282SShri Abhyankar       break;
8922b4ed282SShri Abhyankar     }
8932b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
8942b4ed282SShri Abhyankar     if (snes->domainerror) {
8952b4ed282SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8962b4ed282SShri Abhyankar       PetscFunctionReturn(0);
8972b4ed282SShri Abhyankar     }
898c2fc9fa9SBarry Smith     if (snes->ops->solve != SNESSolve_VISS) {
8997fe79bc4SShri Abhyankar       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
9007fe79bc4SShri Abhyankar     } else {
9012b4ed282SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
9027fe79bc4SShri Abhyankar     }
90362cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
9040661309fSPeter Brune     if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* sufficient reduction */
9050661309fSPeter Brune       if (snes->ls_monitor) {
9060661309fSPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
90762d1f40fSBarry Smith         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Quadratically determined step, lambda=%g\n",(double)lambda);CHKERRQ(ierr);
9080661309fSPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
9092b4ed282SShri Abhyankar       }
9102b4ed282SShri Abhyankar       break;
9112b4ed282SShri Abhyankar     }
9122b4ed282SShri Abhyankar     count++;
9132b4ed282SShri Abhyankar   }
9142b4ed282SShri Abhyankar   theend2:
9152b4ed282SShri Abhyankar   /* Optional user-defined check for line search step validity */
9160661309fSPeter Brune   if (snes->ops->postcheckstep) {
9170661309fSPeter Brune     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
9182b4ed282SShri Abhyankar     if (changed_y) {
9192b4ed282SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
9207fe79bc4SShri Abhyankar       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
9212b4ed282SShri Abhyankar     }
9222b4ed282SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
9232b4ed282SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);
9242b4ed282SShri Abhyankar       if (snes->domainerror) {
9252b4ed282SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
9262b4ed282SShri Abhyankar         PetscFunctionReturn(0);
9272b4ed282SShri Abhyankar       }
928c2fc9fa9SBarry Smith       if (snes->ops->solve != SNESSolve_VISS) {
9297fe79bc4SShri Abhyankar         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
9307fe79bc4SShri Abhyankar       } else {
9317fe79bc4SShri Abhyankar         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
9327fe79bc4SShri Abhyankar       }
9337fe79bc4SShri Abhyankar 
9342b4ed282SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
9352b4ed282SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
93662cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
9372b4ed282SShri Abhyankar     }
9382b4ed282SShri Abhyankar   }
9392b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
9402b4ed282SShri Abhyankar   PetscFunctionReturn(0);
9412b4ed282SShri Abhyankar }
9422b4ed282SShri Abhyankar 
9432b4ed282SShri Abhyankar 
9445559b345SBarry Smith #undef __FUNCT__
9455559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds"
9465559b345SBarry Smith /*@
9472b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
9482b4ed282SShri Abhyankar 
9492b4ed282SShri Abhyankar    Input Parameters:
9502b4ed282SShri Abhyankar .  snes - the SNES context.
9512b4ed282SShri Abhyankar .  xl   - lower bound.
9522b4ed282SShri Abhyankar .  xu   - upper bound.
9532b4ed282SShri Abhyankar 
9542b4ed282SShri Abhyankar    Notes:
9552b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
95601f0b76bSBarry Smith    SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp().
95784c105d7SBarry Smith 
9582bd2b0e6SSatish Balay    Level: advanced
9592bd2b0e6SSatish Balay 
9605559b345SBarry Smith @*/
96169c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
9622b4ed282SShri Abhyankar {
96361589011SJed Brown   PetscErrorCode   ierr,(*f)(SNES,Vec,Vec);
9642b4ed282SShri Abhyankar 
9652b4ed282SShri Abhyankar   PetscFunctionBegin;
9662b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
9672b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
9682b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
96961589011SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
97061589011SJed Brown   if (!f) {ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);}
97161589011SJed Brown   ierr = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr);
97261589011SJed Brown   PetscFunctionReturn(0);
97361589011SJed Brown }
97461589011SJed Brown 
97561589011SJed Brown EXTERN_C_BEGIN
97661589011SJed Brown #undef __FUNCT__
97761589011SJed Brown #define __FUNCT__ "SNESVISetVariableBounds_VI"
97861589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
97961589011SJed Brown {
98061589011SJed Brown   PetscErrorCode    ierr;
98161589011SJed Brown   const PetscScalar *xxl,*xxu;
98261589011SJed Brown   PetscInt          i,n, cnt = 0;
98361589011SJed Brown 
98461589011SJed Brown   PetscFunctionBegin;
985a63bb30eSJed Brown   ierr = SNESGetFunction(snes,&snes->vec_func,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
986a63bb30eSJed Brown   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
9872b4ed282SShri 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);
9882b4ed282SShri 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);
989c2fc9fa9SBarry Smith   ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);
9902176524cSBarry Smith   ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
9912176524cSBarry Smith   ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
992c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
993c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
994c2fc9fa9SBarry Smith   snes->xl = xl;
995c2fc9fa9SBarry Smith   snes->xu = xu;
9966fd67740SBarry Smith   ierr = VecGetLocalSize(xl,&n);CHKERRQ(ierr);
9976fd67740SBarry Smith   ierr = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr);
9986fd67740SBarry Smith   ierr = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr);
9996fd67740SBarry Smith   for (i=0; i<n; i++) {
10006fd67740SBarry Smith     cnt += ((xxl[i] != SNES_VI_NINF) || (xxu[i] != SNES_VI_INF));
10016fd67740SBarry Smith   }
1002c2fc9fa9SBarry Smith   ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
10036fd67740SBarry Smith   ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr);
10046fd67740SBarry Smith   ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr);
10052b4ed282SShri Abhyankar   PetscFunctionReturn(0);
10062b4ed282SShri Abhyankar }
100761589011SJed Brown EXTERN_C_END
100892c02d66SPeter Brune 
100992c02d66SPeter Brune EXTERN_C_BEGIN
101092c02d66SPeter Brune #undef __FUNCT__
101192c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_VI"
101292c02d66SPeter Brune PetscErrorCode  SNESLineSearchSetType_VI(SNES snes, SNESLineSearchType type)
101392c02d66SPeter Brune {
101492c02d66SPeter Brune   PetscErrorCode ierr;
101592c02d66SPeter Brune   PetscFunctionBegin;
101692c02d66SPeter Brune 
101792c02d66SPeter Brune   switch (type) {
101892c02d66SPeter Brune   case SNES_LS_BASIC:
101992c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
102092c02d66SPeter Brune     break;
102192c02d66SPeter Brune   case SNES_LS_BASIC_NONORMS:
102292c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
102392c02d66SPeter Brune     break;
102492c02d66SPeter Brune   case SNES_LS_QUADRATIC:
102592c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
102692c02d66SPeter Brune     break;
102792c02d66SPeter Brune   case SNES_LS_CUBIC:
102892c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
102992c02d66SPeter Brune     break;
103092c02d66SPeter Brune   default:
103192c02d66SPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
103292c02d66SPeter Brune     break;
103392c02d66SPeter Brune   }
103492c02d66SPeter Brune   snes->ls_type = type;
103592c02d66SPeter Brune   PetscFunctionReturn(0);
103692c02d66SPeter Brune }
103792c02d66SPeter Brune EXTERN_C_END
103892c02d66SPeter Brune 
10392b4ed282SShri Abhyankar #undef __FUNCT__
1040c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VI"
1041c2fc9fa9SBarry Smith PetscErrorCode SNESSetFromOptions_VI(SNES snes)
10422b4ed282SShri Abhyankar {
10432b4ed282SShri Abhyankar   PetscErrorCode  ierr;
1044c2fc9fa9SBarry Smith   PetscBool       flg;
10452b4ed282SShri Abhyankar 
10462b4ed282SShri Abhyankar   PetscFunctionBegin;
1047c2fc9fa9SBarry Smith   ierr = PetscOptionsHead("SNES VI options");CHKERRQ(ierr);
1048c2fc9fa9SBarry Smith   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
1049c2fc9fa9SBarry Smith   if (flg) {
1050c2fc9fa9SBarry Smith     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
10512b4ed282SShri Abhyankar   }
1052c2fc9fa9SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1053ed70e96aSJungho Lee   PetscFunctionReturn(0);
1054ed70e96aSJungho Lee }
1055