xref: /petsc/src/snes/impls/vi/vi.c (revision a63bb30eff8b1dff6ddce5531d6230e04f6fd10b)
12b4ed282SShri Abhyankar 
2*c2fc9fa9SBarry 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 {
212176524cSBarry Smith   PetscErrorCode   ierr;
222176524cSBarry Smith 
232176524cSBarry Smith   PetscFunctionBegin;
24*c2fc9fa9SBarry Smith   ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);
25*c2fc9fa9SBarry Smith   snes->ops->computevariablebounds = compute;
262176524cSBarry Smith   PetscFunctionReturn(0);
272176524cSBarry Smith }
282176524cSBarry Smith 
292176524cSBarry Smith 
302176524cSBarry Smith #undef __FUNCT__
31ed70e96aSJungho Lee #define __FUNCT__ "SNESVIComputeInactiveSetIS"
32d0af7cd3SBarry Smith /*
33ed70e96aSJungho Lee    SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables
34d0af7cd3SBarry Smith 
35d0af7cd3SBarry Smith    Input parameter
36d0af7cd3SBarry Smith .  snes - the SNES context
37d0af7cd3SBarry Smith .  X    - the snes solution vector
38d0af7cd3SBarry Smith 
39d0af7cd3SBarry Smith    Output parameter
40d0af7cd3SBarry Smith .  ISact - active set index set
41d0af7cd3SBarry Smith 
42d0af7cd3SBarry Smith  */
43ed70e96aSJungho Lee PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS* inact)
44d0af7cd3SBarry Smith {
45d0af7cd3SBarry Smith   PetscErrorCode    ierr;
46d0af7cd3SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
47d0af7cd3SBarry Smith   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
48d0af7cd3SBarry Smith 
49d0af7cd3SBarry Smith   PetscFunctionBegin;
50d0af7cd3SBarry Smith   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
51d0af7cd3SBarry Smith   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
52d0af7cd3SBarry Smith   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
53d0af7cd3SBarry Smith   ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr);
54d0af7cd3SBarry Smith   ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr);
55d0af7cd3SBarry Smith   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
56d0af7cd3SBarry Smith   /* Compute inactive set size */
57d0af7cd3SBarry Smith   for (i=0; i < nlocal;i++) {
58d0af7cd3SBarry 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++;
59d0af7cd3SBarry Smith   }
60d0af7cd3SBarry Smith 
61d0af7cd3SBarry Smith   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
62d0af7cd3SBarry Smith 
63d0af7cd3SBarry Smith   /* Set inactive set indices */
64d0af7cd3SBarry Smith   for(i=0; i < nlocal; i++) {
65d0af7cd3SBarry 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;
66d0af7cd3SBarry Smith   }
67d0af7cd3SBarry Smith 
68d0af7cd3SBarry Smith    /* Create inactive set IS */
69d0af7cd3SBarry Smith   ierr = ISCreateGeneral(((PetscObject)upper)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr);
70d0af7cd3SBarry Smith 
71d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
72d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr);
73d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr);
74d0af7cd3SBarry Smith   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
75d0af7cd3SBarry Smith   PetscFunctionReturn(0);
76d0af7cd3SBarry Smith }
77d0af7cd3SBarry Smith 
783c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/
792b4ed282SShri Abhyankar 
809308c137SBarry Smith #undef __FUNCT__
819308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI"
827087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
839308c137SBarry Smith {
849308c137SBarry Smith   PetscErrorCode    ierr;
85649052a6SBarry Smith   PetscViewer        viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
869308c137SBarry Smith   const PetscScalar  *x,*xl,*xu,*f;
876fd67740SBarry Smith   PetscInt           i,n,act[2] = {0,0},fact[2],N;
886a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
896a9e2d46SJungho Lee   PetscInt           act_bound[2] = {0,0},fact_bound[2];
909308c137SBarry Smith   PetscReal          rnorm,fnorm;
919308c137SBarry Smith 
929308c137SBarry Smith   PetscFunctionBegin;
939308c137SBarry Smith   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
946fd67740SBarry Smith   ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr);
95*c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
96*c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
979308c137SBarry Smith   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
983f731af5SBarry Smith   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
999308c137SBarry Smith 
1009308c137SBarry Smith   rnorm = 0.0;
1019308c137SBarry Smith   for (i=0; i<n; i++) {
10210f5ae6bSBarry 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]);
103e400df7aSBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++;
104e400df7aSBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++;
105e400df7aSBarry Smith     else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Can never get here");
1069308c137SBarry Smith   }
1076a9e2d46SJungho Lee 
1086a9e2d46SJungho Lee   for (i=0; i<n; i++) {
1096a9e2d46SJungho Lee     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++;
1106a9e2d46SJungho Lee     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++;
1116a9e2d46SJungho Lee   }
1123f731af5SBarry Smith   ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
113*c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
114*c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
1159308c137SBarry Smith   ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
116d9822059SBarry Smith   ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
11721a4c9b1SBarry Smith   ierr = MPI_Allreduce(act,fact,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
1186a9e2d46SJungho Lee   ierr = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
119f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1206fd67740SBarry Smith 
121649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
122*c2fc9fa9SBarry 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);
1236a9e2d46SJungho Lee 
124649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1259308c137SBarry Smith   PetscFunctionReturn(0);
1269308c137SBarry Smith }
1279308c137SBarry Smith 
1282b4ed282SShri Abhyankar /*
1292b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1302b4ed282SShri 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
1312b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1322b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1332b4ed282SShri Abhyankar */
1342b4ed282SShri Abhyankar #undef __FUNCT__
1352b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private"
136ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
1372b4ed282SShri Abhyankar {
1382b4ed282SShri Abhyankar   PetscReal      a1;
1392b4ed282SShri Abhyankar   PetscErrorCode ierr;
140ace3abfcSBarry Smith   PetscBool     hastranspose;
1412b4ed282SShri Abhyankar 
1422b4ed282SShri Abhyankar   PetscFunctionBegin;
1432b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1442b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1452b4ed282SShri Abhyankar   if (hastranspose) {
1462b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1472b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
1482b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
1494839bfe8SBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
1502b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
1512b4ed282SShri Abhyankar   } else {
1522b4ed282SShri Abhyankar     Vec         work;
1532b4ed282SShri Abhyankar     PetscScalar result;
1542b4ed282SShri Abhyankar     PetscReal   wnorm;
1552b4ed282SShri Abhyankar 
1562b4ed282SShri Abhyankar     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
1572b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
1582b4ed282SShri Abhyankar     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
1592b4ed282SShri Abhyankar     ierr = MatMult(A,W,work);CHKERRQ(ierr);
1602b4ed282SShri Abhyankar     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
1616bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
1622b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
1634839bfe8SBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
1642b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1652b4ed282SShri Abhyankar   }
1662b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1672b4ed282SShri Abhyankar }
1682b4ed282SShri Abhyankar 
1692b4ed282SShri Abhyankar /*
1702b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
1712b4ed282SShri Abhyankar */
1722b4ed282SShri Abhyankar #undef __FUNCT__
1732b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private"
1742b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
1752b4ed282SShri Abhyankar {
1762b4ed282SShri Abhyankar   PetscReal      a1,a2;
1772b4ed282SShri Abhyankar   PetscErrorCode ierr;
178ace3abfcSBarry Smith   PetscBool     hastranspose;
1792b4ed282SShri Abhyankar 
1802b4ed282SShri Abhyankar   PetscFunctionBegin;
1812b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1822b4ed282SShri Abhyankar   if (hastranspose) {
1832b4ed282SShri Abhyankar     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
1842b4ed282SShri Abhyankar     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
1852b4ed282SShri Abhyankar 
1862b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1872b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
1882b4ed282SShri Abhyankar     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
1892b4ed282SShri Abhyankar     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
1902b4ed282SShri Abhyankar     if (a1 != 0.0) {
1914839bfe8SBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
1922b4ed282SShri Abhyankar     }
1932b4ed282SShri Abhyankar   }
1942b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1952b4ed282SShri Abhyankar }
1962b4ed282SShri Abhyankar 
1972b4ed282SShri Abhyankar /*
1982b4ed282SShri Abhyankar   SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm.
1992b4ed282SShri Abhyankar 
2002b4ed282SShri Abhyankar   Notes:
2012b4ed282SShri Abhyankar   The convergence criterion currently implemented is
2022b4ed282SShri Abhyankar   merit < abstol
2032b4ed282SShri Abhyankar   merit < rtol*merit_initial
2042b4ed282SShri Abhyankar */
2052b4ed282SShri Abhyankar #undef __FUNCT__
2062b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI"
2077fe79bc4SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
2082b4ed282SShri Abhyankar {
2092b4ed282SShri Abhyankar   PetscErrorCode ierr;
2102b4ed282SShri Abhyankar 
2112b4ed282SShri Abhyankar   PetscFunctionBegin;
2122b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2132b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
2142b4ed282SShri Abhyankar 
2152b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
2162b4ed282SShri Abhyankar 
2172b4ed282SShri Abhyankar   if (!it) {
2182b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
2197fe79bc4SShri Abhyankar     snes->ttol = fnorm*snes->rtol;
2202b4ed282SShri Abhyankar   }
2217fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
2222b4ed282SShri Abhyankar     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
2232b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
2247fe79bc4SShri Abhyankar   } else if (fnorm < snes->abstol) {
2254839bfe8SBarry Smith     ierr = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
2262b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
2272b4ed282SShri Abhyankar   } else if (snes->nfuncs >= snes->max_funcs) {
2282b4ed282SShri Abhyankar     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
2292b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
2302b4ed282SShri Abhyankar   }
2312b4ed282SShri Abhyankar 
2322b4ed282SShri Abhyankar   if (it && !*reason) {
2337fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2344839bfe8SBarry Smith       ierr = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
2352b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2362b4ed282SShri Abhyankar     }
2372b4ed282SShri Abhyankar   }
2382b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2392b4ed282SShri Abhyankar }
2402b4ed282SShri Abhyankar 
241ee657d29SShri Abhyankar 
242c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
243c1a5e521SShri Abhyankar /*
244c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
245c1a5e521SShri Abhyankar 
246c1a5e521SShri Abhyankar    Input Parameters:
247c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
248c1a5e521SShri Abhyankar 
249c1a5e521SShri Abhyankar    Output Parameters:
250c1a5e521SShri Abhyankar .  X - Bound projected X
251c1a5e521SShri Abhyankar 
252c1a5e521SShri Abhyankar */
253c1a5e521SShri Abhyankar 
254c1a5e521SShri Abhyankar #undef __FUNCT__
255c1a5e521SShri Abhyankar #define __FUNCT__ "SNESVIProjectOntoBounds"
256c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
257c1a5e521SShri Abhyankar {
258c1a5e521SShri Abhyankar   PetscErrorCode    ierr;
259c1a5e521SShri Abhyankar   const PetscScalar *xl,*xu;
260c1a5e521SShri Abhyankar   PetscScalar       *x;
261c1a5e521SShri Abhyankar   PetscInt          i,n;
262c1a5e521SShri Abhyankar 
263c1a5e521SShri Abhyankar   PetscFunctionBegin;
264c1a5e521SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
265c1a5e521SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
266*c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
267*c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
268c1a5e521SShri Abhyankar 
269c1a5e521SShri Abhyankar   for(i = 0;i<n;i++) {
27010f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
27110f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
272c1a5e521SShri Abhyankar   }
273c1a5e521SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
274*c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
275*c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
276c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
277c1a5e521SShri Abhyankar }
278c1a5e521SShri Abhyankar 
27969c03620SShri Abhyankar 
28069c03620SShri Abhyankar #undef __FUNCT__
281693ddf92SShri Abhyankar #define __FUNCT__ "SNESVIGetActiveSetIS"
282693ddf92SShri Abhyankar /*
283693ddf92SShri Abhyankar    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
284693ddf92SShri Abhyankar 
285693ddf92SShri Abhyankar    Input parameter
286693ddf92SShri Abhyankar .  snes - the SNES context
287693ddf92SShri Abhyankar .  X    - the snes solution vector
288693ddf92SShri Abhyankar .  F    - the nonlinear function vector
289693ddf92SShri Abhyankar 
290693ddf92SShri Abhyankar    Output parameter
291693ddf92SShri Abhyankar .  ISact - active set index set
292693ddf92SShri Abhyankar  */
293693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS* ISact)
294d950fb63SShri Abhyankar {
295d950fb63SShri Abhyankar   PetscErrorCode   ierr;
296*c2fc9fa9SBarry Smith   Vec               Xl=snes->xl,Xu=snes->xu;
297693ddf92SShri Abhyankar   const PetscScalar *x,*f,*xl,*xu;
298693ddf92SShri Abhyankar   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
299d950fb63SShri Abhyankar 
300d950fb63SShri Abhyankar   PetscFunctionBegin;
301d950fb63SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
302d950fb63SShri Abhyankar   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
303fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
304fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
305fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
306fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
307693ddf92SShri Abhyankar   /* Compute active set size */
308d950fb63SShri Abhyankar   for (i=0; i < nlocal;i++) {
30910f5ae6bSBarry 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++;
310d950fb63SShri Abhyankar   }
311693ddf92SShri Abhyankar 
312d950fb63SShri Abhyankar   ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr);
313d950fb63SShri Abhyankar 
314693ddf92SShri Abhyankar   /* Set active set indices */
315d950fb63SShri Abhyankar   for(i=0; i < nlocal; i++) {
31610f5ae6bSBarry 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;
317d950fb63SShri Abhyankar   }
318d950fb63SShri Abhyankar 
319693ddf92SShri Abhyankar    /* Create active set IS */
32062298a1eSBarry Smith   ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
321d950fb63SShri Abhyankar 
322fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
323fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
324fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
325fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
326d950fb63SShri Abhyankar   PetscFunctionReturn(0);
327d950fb63SShri Abhyankar }
328d950fb63SShri Abhyankar 
329693ddf92SShri Abhyankar #undef __FUNCT__
330693ddf92SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS"
331693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS* ISact,IS* ISinact)
332693ddf92SShri Abhyankar {
333693ddf92SShri Abhyankar   PetscErrorCode     ierr;
334693ddf92SShri Abhyankar 
335693ddf92SShri Abhyankar   PetscFunctionBegin;
336693ddf92SShri Abhyankar   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
337693ddf92SShri Abhyankar   ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr);
338693ddf92SShri Abhyankar   PetscFunctionReturn(0);
339693ddf92SShri Abhyankar }
340693ddf92SShri Abhyankar 
341fe835674SShri Abhyankar #undef __FUNCT__
342fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
34310f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X,PetscReal *fnorm)
344fe835674SShri Abhyankar {
345fe835674SShri Abhyankar   PetscErrorCode    ierr;
346fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
347fe835674SShri Abhyankar   PetscInt          i,n;
348fe835674SShri Abhyankar   PetscReal         rnorm;
349fe835674SShri Abhyankar 
350fe835674SShri Abhyankar   PetscFunctionBegin;
351fe835674SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
352*c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
353*c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
354fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
355fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
356fe835674SShri Abhyankar   rnorm = 0.0;
357fe835674SShri Abhyankar   for (i=0; i<n; i++) {
35810f5ae6bSBarry 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]);
359fe835674SShri Abhyankar   }
360fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
361*c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
362*c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
363fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
364d9822059SBarry Smith   ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
36562d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
366fe835674SShri Abhyankar   PetscFunctionReturn(0);
367fe835674SShri Abhyankar }
368fe835674SShri Abhyankar 
3692f63a38cSShri Abhyankar 
3702b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3712b4ed282SShri Abhyankar /*
372*c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3732b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3742b4ed282SShri Abhyankar 
3752b4ed282SShri Abhyankar    Input Parameter:
3762b4ed282SShri Abhyankar .  snes - the SNES context
3772b4ed282SShri Abhyankar .  x - the solution vector
3782b4ed282SShri Abhyankar 
3792b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
3802b4ed282SShri Abhyankar 
3812b4ed282SShri Abhyankar    Notes:
3822b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
3832b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
3842b4ed282SShri Abhyankar    the call to SNESSolve().
3852b4ed282SShri Abhyankar  */
3862b4ed282SShri Abhyankar #undef __FUNCT__
3872b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
3882b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
3892b4ed282SShri Abhyankar {
3902b4ed282SShri Abhyankar   PetscErrorCode ierr;
3912b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
3922b4ed282SShri Abhyankar 
3932b4ed282SShri Abhyankar   PetscFunctionBegin;
39458c9b817SLisandro Dalcin 
39558c9b817SLisandro Dalcin   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
3962b4ed282SShri Abhyankar 
397*c2fc9fa9SBarry Smith   if (snes->ops->computevariablebounds) {
398*c2fc9fa9SBarry Smith     if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
399*c2fc9fa9SBarry Smith     if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
400*c2fc9fa9SBarry Smith     ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
401*c2fc9fa9SBarry Smith   } else if (!snes->xl && !snes->xu) {
4022176524cSBarry Smith     /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
403*c2fc9fa9SBarry Smith     ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr);
404*c2fc9fa9SBarry Smith     ierr = VecSet(snes->xl,SNES_VI_NINF);CHKERRQ(ierr);
405*c2fc9fa9SBarry Smith     ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr);
406*c2fc9fa9SBarry Smith     ierr = VecSet(snes->xu,SNES_VI_INF);CHKERRQ(ierr);
4072b4ed282SShri Abhyankar   } else {
4082b4ed282SShri Abhyankar     /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
4092b4ed282SShri Abhyankar     ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
410*c2fc9fa9SBarry Smith     ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr);
411*c2fc9fa9SBarry Smith     ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr);
4122b4ed282SShri 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]))
4132b4ed282SShri 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.");
4142b4ed282SShri Abhyankar   }
415abcba341SShri Abhyankar 
4162b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4172b4ed282SShri Abhyankar }
4182b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
4192176524cSBarry Smith #undef __FUNCT__
4202176524cSBarry Smith #define __FUNCT__ "SNESReset_VI"
4212176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes)
4222176524cSBarry Smith {
4232176524cSBarry Smith   PetscErrorCode ierr;
4242176524cSBarry Smith 
4252176524cSBarry Smith   PetscFunctionBegin;
426*c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
427*c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
4282176524cSBarry Smith   PetscFunctionReturn(0);
4292176524cSBarry Smith }
4302176524cSBarry Smith 
4312b4ed282SShri Abhyankar /*
4322b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
4332b4ed282SShri Abhyankar    with SNESCreate_VI().
4342b4ed282SShri Abhyankar 
4352b4ed282SShri Abhyankar    Input Parameter:
4362b4ed282SShri Abhyankar .  snes - the SNES context
4372b4ed282SShri Abhyankar 
4382b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4392b4ed282SShri Abhyankar  */
4402b4ed282SShri Abhyankar #undef __FUNCT__
4412b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI"
4422b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
4432b4ed282SShri Abhyankar {
4442b4ed282SShri Abhyankar   PetscErrorCode ierr;
4452b4ed282SShri Abhyankar 
4462b4ed282SShri Abhyankar   PetscFunctionBegin;
4472b4ed282SShri Abhyankar   ierr = PetscFree(snes->data);CHKERRQ(ierr);
4482b4ed282SShri Abhyankar 
4492b4ed282SShri Abhyankar   /* clear composed functions */
4502b4ed282SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
451c92abb8aSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
4522b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4532b4ed282SShri Abhyankar }
4547fe79bc4SShri Abhyankar 
4552b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
456*c2fc9fa9SBarry Smith 
457*c2fc9fa9SBarry Smith /*
458*c2fc9fa9SBarry Smith 
459*c2fc9fa9SBarry Smith       These line searches are common for all the VI solvers
460*c2fc9fa9SBarry Smith */
461*c2fc9fa9SBarry Smith extern PetscErrorCode SNESSolve_VISS(SNES);
462*c2fc9fa9SBarry Smith 
4632b4ed282SShri Abhyankar #undef __FUNCT__
4642b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI"
4652b4ed282SShri Abhyankar 
4662b4ed282SShri Abhyankar /*
4677fe79bc4SShri Abhyankar   This routine does not actually do a line search but it takes a full newton
4687fe79bc4SShri Abhyankar   step while ensuring that the new iterates remain within the constraints.
4692b4ed282SShri Abhyankar 
4702b4ed282SShri Abhyankar */
4711e633543SBarry 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)
4722b4ed282SShri Abhyankar {
4732b4ed282SShri Abhyankar   PetscErrorCode ierr;
474ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
4752b4ed282SShri Abhyankar 
4762b4ed282SShri Abhyankar   PetscFunctionBegin;
4772b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
4782b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4792b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
4802b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
481c1a5e521SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
4820661309fSPeter Brune   if (snes->ops->postcheckstep) {
4830661309fSPeter Brune    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
484c1a5e521SShri Abhyankar   }
485c1a5e521SShri Abhyankar   if (changed_y) {
486c1a5e521SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4877fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
488c1a5e521SShri Abhyankar   }
489c1a5e521SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
490c1a5e521SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
491c1a5e521SShri Abhyankar   if (!snes->domainerror) {
492*c2fc9fa9SBarry Smith     if (snes->ops->solve != SNESSolve_VISS) {
4937fe79bc4SShri Abhyankar        ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
4947fe79bc4SShri Abhyankar     } else {
495c1a5e521SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
4967fe79bc4SShri Abhyankar     }
49762cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
498c1a5e521SShri Abhyankar   }
4990661309fSPeter Brune   if (snes->ls_monitor) {
5000661309fSPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
5010661309fSPeter Brune     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
5020661309fSPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
503c1a5e521SShri Abhyankar   }
504c1a5e521SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
505c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
506c1a5e521SShri Abhyankar }
507c1a5e521SShri Abhyankar 
508c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
509c1a5e521SShri Abhyankar #undef __FUNCT__
5102b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI"
5112b4ed282SShri Abhyankar 
5122b4ed282SShri Abhyankar /*
5132b4ed282SShri Abhyankar   This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c
5142b4ed282SShri Abhyankar */
5151e633543SBarry 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)
5162b4ed282SShri Abhyankar {
5172b4ed282SShri Abhyankar   PetscErrorCode ierr;
518ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5192b4ed282SShri Abhyankar 
5202b4ed282SShri Abhyankar   PetscFunctionBegin;
5212b4ed282SShri Abhyankar   *flag = PETSC_TRUE;
5222b4ed282SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5232b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
5247fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
5250661309fSPeter Brune   if (snes->ops->postcheckstep) {
5260661309fSPeter Brune    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
5272b4ed282SShri Abhyankar   }
5282b4ed282SShri Abhyankar   if (changed_y) {
5292b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
5307fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
5312b4ed282SShri Abhyankar   }
5322b4ed282SShri Abhyankar 
5332b4ed282SShri Abhyankar   /* don't evaluate function the last time through */
5342b4ed282SShri Abhyankar   if (snes->iter < snes->max_its-1) {
5352b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
5362b4ed282SShri Abhyankar   }
5372b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5382b4ed282SShri Abhyankar   PetscFunctionReturn(0);
5392b4ed282SShri Abhyankar }
5407fe79bc4SShri Abhyankar 
5412b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
542*c2fc9fa9SBarry Smith 
5432b4ed282SShri Abhyankar #undef __FUNCT__
5442b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI"
5452b4ed282SShri Abhyankar /*
5467fe79bc4SShri Abhyankar   This routine implements a cubic line search while doing a projection on the variable bounds
5472b4ed282SShri Abhyankar */
5481e633543SBarry 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)
5492b4ed282SShri Abhyankar {
550fe835674SShri Abhyankar   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
551fe835674SShri Abhyankar   PetscReal      minlambda,lambda,lambdatemp;
552fe835674SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
553fe835674SShri Abhyankar   PetscScalar    cinitslope;
554fe835674SShri Abhyankar #endif
555fe835674SShri Abhyankar   PetscErrorCode ierr;
556fe835674SShri Abhyankar   PetscInt       count;
557fe835674SShri Abhyankar   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
558fe835674SShri Abhyankar   MPI_Comm       comm;
559fe835674SShri Abhyankar 
560fe835674SShri Abhyankar   PetscFunctionBegin;
561fe835674SShri Abhyankar   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
562fe835674SShri Abhyankar   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
563fe835674SShri Abhyankar   *flag   = PETSC_TRUE;
564fe835674SShri Abhyankar 
565fe835674SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
566fe835674SShri Abhyankar   if (*ynorm == 0.0) {
5670661309fSPeter Brune     if (snes->ls_monitor) {
5680661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
5690661309fSPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
5700661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
571fe835674SShri Abhyankar     }
572fe835674SShri Abhyankar     *gnorm = fnorm;
573fe835674SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
574fe835674SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
575fe835674SShri Abhyankar     *flag  = PETSC_FALSE;
576fe835674SShri Abhyankar     goto theend1;
577fe835674SShri Abhyankar   }
5780661309fSPeter Brune   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
5790661309fSPeter Brune     if (snes->ls_monitor) {
5800661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
58162d1f40fSBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Scaling step by %g old ynorm %g\n",(double)(snes->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr);
5820661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
583fe835674SShri Abhyankar     }
5840661309fSPeter Brune     ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
5850661309fSPeter Brune     *ynorm = snes->maxstep;
586fe835674SShri Abhyankar   }
587fe835674SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5880661309fSPeter Brune   minlambda = snes->steptol/rellength;
589fe835674SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
590fe835674SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
591fe835674SShri Abhyankar   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
592fe835674SShri Abhyankar   initslope = PetscRealPart(cinitslope);
593fe835674SShri Abhyankar #else
594fe835674SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
595fe835674SShri Abhyankar #endif
596fe835674SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
597fe835674SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
598fe835674SShri Abhyankar 
599fe835674SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
600fe835674SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
601fe835674SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
602fe835674SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
603fe835674SShri Abhyankar     *flag = PETSC_FALSE;
604fe835674SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
605fe835674SShri Abhyankar     goto theend1;
606fe835674SShri Abhyankar   }
607fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
608*c2fc9fa9SBarry Smith   if (snes->ops->solve != SNESSolve_VISS) {
609fe835674SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
6107fe79bc4SShri Abhyankar   } else {
6117fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
6127fe79bc4SShri Abhyankar   }
613fe835674SShri Abhyankar   if (snes->domainerror) {
614fe835674SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
615fe835674SShri Abhyankar     PetscFunctionReturn(0);
616fe835674SShri Abhyankar   }
61762cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6180661309fSPeter 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);
6190661309fSPeter Brune   if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* Sufficient reduction */
6200661309fSPeter Brune     if (snes->ls_monitor) {
6210661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
62262d1f40fSBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr);
6230661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
624fe835674SShri Abhyankar     }
625fe835674SShri Abhyankar     goto theend1;
626fe835674SShri Abhyankar   }
627fe835674SShri Abhyankar 
628fe835674SShri Abhyankar   /* Fit points with quadratic */
629fe835674SShri Abhyankar   lambda     = 1.0;
630fe835674SShri Abhyankar   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
631fe835674SShri Abhyankar   lambdaprev = lambda;
632fe835674SShri Abhyankar   gnormprev  = *gnorm;
633fe835674SShri Abhyankar   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
634fe835674SShri Abhyankar   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
635fe835674SShri Abhyankar   else                         lambda = lambdatemp;
636fe835674SShri Abhyankar 
637fe835674SShri Abhyankar   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
638fe835674SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
639fe835674SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
640fe835674SShri Abhyankar     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
641fe835674SShri Abhyankar     *flag = PETSC_FALSE;
642fe835674SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
643fe835674SShri Abhyankar     goto theend1;
644fe835674SShri Abhyankar   }
645fe835674SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
646*c2fc9fa9SBarry Smith   if (snes->ops->solve != SNESSolve_VISS) {
647fe835674SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
6487fe79bc4SShri Abhyankar   } else {
6497fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
6507fe79bc4SShri Abhyankar   }
651fe835674SShri Abhyankar   if (snes->domainerror) {
652fe835674SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
653fe835674SShri Abhyankar     PetscFunctionReturn(0);
654fe835674SShri Abhyankar   }
65562cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6560661309fSPeter Brune   if (snes->ls_monitor) {
6570661309fSPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
65862d1f40fSBarry Smith     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: gnorm after quadratic fit %g\n",(double)(*gnorm));CHKERRQ(ierr);
6590661309fSPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
660fe835674SShri Abhyankar   }
6610661309fSPeter Brune   if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* sufficient reduction */
6620661309fSPeter Brune     if (snes->ls_monitor) {
6630661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
6644839bfe8SBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
6650661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
666fe835674SShri Abhyankar     }
667fe835674SShri Abhyankar     goto theend1;
668fe835674SShri Abhyankar   }
669fe835674SShri Abhyankar 
670fe835674SShri Abhyankar   /* Fit points with cubic */
671fe835674SShri Abhyankar   count = 1;
672fe835674SShri Abhyankar   while (PETSC_TRUE) {
673fe835674SShri Abhyankar     if (lambda <= minlambda) {
6740661309fSPeter Brune       if (snes->ls_monitor) {
6750661309fSPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
6760661309fSPeter Brune  	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
67762d1f40fSBarry 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);
6780661309fSPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
679fe835674SShri Abhyankar       }
680fe835674SShri Abhyankar       *flag = PETSC_FALSE;
681fe835674SShri Abhyankar       break;
682fe835674SShri Abhyankar     }
683fe835674SShri Abhyankar     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
684fe835674SShri Abhyankar     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
685fe835674SShri Abhyankar     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
686fe835674SShri Abhyankar     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
687fe835674SShri Abhyankar     d  = b*b - 3*a*initslope;
688fe835674SShri Abhyankar     if (d < 0.0) d = 0.0;
689fe835674SShri Abhyankar     if (a == 0.0) {
690fe835674SShri Abhyankar       lambdatemp = -initslope/(2.0*b);
691fe835674SShri Abhyankar     } else {
69262d1f40fSBarry Smith       lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
693fe835674SShri Abhyankar     }
694fe835674SShri Abhyankar     lambdaprev = lambda;
695fe835674SShri Abhyankar     gnormprev  = *gnorm;
696fe835674SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
697fe835674SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
698fe835674SShri Abhyankar     else                         lambda     = lambdatemp;
699fe835674SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
700fe835674SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
701fe835674SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
702fe835674SShri Abhyankar       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
7034839bfe8SBarry 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);
704fe835674SShri Abhyankar       *flag = PETSC_FALSE;
705fe835674SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
706fe835674SShri Abhyankar       break;
707fe835674SShri Abhyankar     }
708fe835674SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
709*c2fc9fa9SBarry Smith     if (snes->ops->solve != SNESSolve_VISS) {
710fe835674SShri Abhyankar       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
7117fe79bc4SShri Abhyankar     } else {
7127fe79bc4SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7137fe79bc4SShri Abhyankar     }
714fe835674SShri Abhyankar     if (snes->domainerror) {
715fe835674SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
716fe835674SShri Abhyankar       PetscFunctionReturn(0);
717fe835674SShri Abhyankar     }
71862cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7190661309fSPeter Brune     if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* is reduction enough? */
7200661309fSPeter Brune       if (snes->ls_monitor) {
72162d1f40fSBarry Smith 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %g lambda=%18.16e\n",(double)(*gnorm),(double)lambda);CHKERRQ(ierr);
722fe835674SShri Abhyankar       }
723fe835674SShri Abhyankar       break;
724fe835674SShri Abhyankar     } else {
7250661309fSPeter Brune       if (snes->ls_monitor) {
72662d1f40fSBarry 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);
727fe835674SShri Abhyankar       }
728fe835674SShri Abhyankar     }
729fe835674SShri Abhyankar     count++;
730fe835674SShri Abhyankar   }
731fe835674SShri Abhyankar   theend1:
732fe835674SShri Abhyankar   /* Optional user-defined check for line search step validity */
7330661309fSPeter Brune   if (snes->ops->postcheckstep && *flag) {
7340661309fSPeter Brune     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
735fe835674SShri Abhyankar     if (changed_y) {
736fe835674SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
737fe835674SShri Abhyankar       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
738fe835674SShri Abhyankar     }
739fe835674SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
740fe835674SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
741*c2fc9fa9SBarry Smith       if (snes->ops->solve != SNESSolve_VISS) {
742fe835674SShri Abhyankar         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
7437fe79bc4SShri Abhyankar       } else {
7447fe79bc4SShri Abhyankar         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7457fe79bc4SShri Abhyankar       }
746fe835674SShri Abhyankar       if (snes->domainerror) {
747fe835674SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
748fe835674SShri Abhyankar         PetscFunctionReturn(0);
749fe835674SShri Abhyankar       }
75062cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
751fe835674SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
752fe835674SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
753fe835674SShri Abhyankar     }
754fe835674SShri Abhyankar   }
755fe835674SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
756fe835674SShri Abhyankar   PetscFunctionReturn(0);
757fe835674SShri Abhyankar }
758fe835674SShri Abhyankar 
7592b4ed282SShri Abhyankar #undef __FUNCT__
7602b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI"
7612b4ed282SShri Abhyankar /*
7627fe79bc4SShri Abhyankar   This routine does a quadratic line search while keeping the iterates within the variable bounds
7632b4ed282SShri Abhyankar */
7641e633543SBarry 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)
7652b4ed282SShri Abhyankar {
7662b4ed282SShri Abhyankar   /*
7672b4ed282SShri Abhyankar      Note that for line search purposes we work with with the related
7682b4ed282SShri Abhyankar      minimization problem:
7692b4ed282SShri Abhyankar         min  z(x):  R^n -> R,
7702b4ed282SShri Abhyankar      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
7712b4ed282SShri Abhyankar    */
7722b4ed282SShri Abhyankar   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
7732b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
7742b4ed282SShri Abhyankar   PetscScalar    cinitslope;
7752b4ed282SShri Abhyankar #endif
7762b4ed282SShri Abhyankar   PetscErrorCode ierr;
7772b4ed282SShri Abhyankar   PetscInt       count;
778ace3abfcSBarry Smith   PetscBool     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7792b4ed282SShri Abhyankar 
7802b4ed282SShri Abhyankar   PetscFunctionBegin;
7812b4ed282SShri Abhyankar   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7822b4ed282SShri Abhyankar   *flag   = PETSC_TRUE;
7832b4ed282SShri Abhyankar 
7842b4ed282SShri Abhyankar   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
7852b4ed282SShri Abhyankar   if (*ynorm == 0.0) {
7860661309fSPeter Brune     if (snes->ls_monitor) {
7870661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
7880661309fSPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
7890661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
7902b4ed282SShri Abhyankar     }
7912b4ed282SShri Abhyankar     *gnorm = fnorm;
7922b4ed282SShri Abhyankar     ierr   = VecCopy(x,w);CHKERRQ(ierr);
7932b4ed282SShri Abhyankar     ierr   = VecCopy(f,g);CHKERRQ(ierr);
7942b4ed282SShri Abhyankar     *flag  = PETSC_FALSE;
7952b4ed282SShri Abhyankar     goto theend2;
7962b4ed282SShri Abhyankar   }
7970661309fSPeter Brune   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
7980661309fSPeter Brune     ierr   = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
7990661309fSPeter Brune     *ynorm = snes->maxstep;
8002b4ed282SShri Abhyankar   }
8012b4ed282SShri Abhyankar   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
8020661309fSPeter Brune   minlambda = snes->steptol/rellength;
8037fe79bc4SShri Abhyankar   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
8042b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
8057fe79bc4SShri Abhyankar   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
8062b4ed282SShri Abhyankar   initslope = PetscRealPart(cinitslope);
8072b4ed282SShri Abhyankar #else
8087fe79bc4SShri Abhyankar   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
8092b4ed282SShri Abhyankar #endif
8102b4ed282SShri Abhyankar   if (initslope > 0.0)  initslope = -initslope;
8112b4ed282SShri Abhyankar   if (initslope == 0.0) initslope = -1.0;
8124839bfe8SBarry Smith   ierr = PetscInfo1(snes,"Initslope %g \n",(double)initslope);CHKERRQ(ierr);
8132b4ed282SShri Abhyankar 
8142b4ed282SShri Abhyankar   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8157fe79bc4SShri Abhyankar   ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
8162b4ed282SShri Abhyankar   if (snes->nfuncs >= snes->max_funcs) {
8172b4ed282SShri Abhyankar     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
8182b4ed282SShri Abhyankar     *flag = PETSC_FALSE;
8192b4ed282SShri Abhyankar     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
8202b4ed282SShri Abhyankar     goto theend2;
8212b4ed282SShri Abhyankar   }
8222b4ed282SShri Abhyankar   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
823*c2fc9fa9SBarry Smith   if (snes->ops->solve != SNESSolve_VISS) {
8247fe79bc4SShri Abhyankar     ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
8257fe79bc4SShri Abhyankar   } else {
8267fe79bc4SShri Abhyankar     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
8277fe79bc4SShri Abhyankar   }
8282b4ed282SShri Abhyankar   if (snes->domainerror) {
8292b4ed282SShri Abhyankar     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8302b4ed282SShri Abhyankar     PetscFunctionReturn(0);
8312b4ed282SShri Abhyankar   }
83262cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8330661309fSPeter Brune   if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* Sufficient reduction */
8340661309fSPeter Brune     if (snes->ls_monitor) {
8350661309fSPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
83662d1f40fSBarry Smith       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr);
8370661309fSPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8382b4ed282SShri Abhyankar     }
8392b4ed282SShri Abhyankar     goto theend2;
8402b4ed282SShri Abhyankar   }
8412b4ed282SShri Abhyankar 
8422b4ed282SShri Abhyankar   /* Fit points with quadratic */
8432b4ed282SShri Abhyankar   lambda = 1.0;
8442b4ed282SShri Abhyankar   count = 1;
8452b4ed282SShri Abhyankar   while (PETSC_TRUE) {
8462b4ed282SShri Abhyankar     if (lambda <= minlambda) { /* bad luck; use full step */
8470661309fSPeter Brune       if (snes->ls_monitor) {
8480661309fSPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8490661309fSPeter Brune         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
85062d1f40fSBarry 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);
8510661309fSPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8522b4ed282SShri Abhyankar       }
8532b4ed282SShri Abhyankar       ierr = VecCopy(x,w);CHKERRQ(ierr);
8542b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
8552b4ed282SShri Abhyankar       break;
8562b4ed282SShri Abhyankar     }
8572b4ed282SShri Abhyankar     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
8582b4ed282SShri Abhyankar     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
8592b4ed282SShri Abhyankar     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
8602b4ed282SShri Abhyankar     else                         lambda     = lambdatemp;
8612b4ed282SShri Abhyankar 
8622b4ed282SShri Abhyankar     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
8637fe79bc4SShri Abhyankar     ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
8642b4ed282SShri Abhyankar     if (snes->nfuncs >= snes->max_funcs) {
8652b4ed282SShri Abhyankar       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
8664839bfe8SBarry 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);
8672b4ed282SShri Abhyankar       *flag = PETSC_FALSE;
8682b4ed282SShri Abhyankar       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
8692b4ed282SShri Abhyankar       break;
8702b4ed282SShri Abhyankar     }
8712b4ed282SShri Abhyankar     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
8722b4ed282SShri Abhyankar     if (snes->domainerror) {
8732b4ed282SShri Abhyankar       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8742b4ed282SShri Abhyankar       PetscFunctionReturn(0);
8752b4ed282SShri Abhyankar     }
876*c2fc9fa9SBarry Smith     if (snes->ops->solve != SNESSolve_VISS) {
8777fe79bc4SShri Abhyankar       ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
8787fe79bc4SShri Abhyankar     } else {
8792b4ed282SShri Abhyankar       ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
8807fe79bc4SShri Abhyankar     }
88162cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8820661309fSPeter Brune     if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* sufficient reduction */
8830661309fSPeter Brune       if (snes->ls_monitor) {
8840661309fSPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
88562d1f40fSBarry Smith         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Quadratically determined step, lambda=%g\n",(double)lambda);CHKERRQ(ierr);
8860661309fSPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8872b4ed282SShri Abhyankar       }
8882b4ed282SShri Abhyankar       break;
8892b4ed282SShri Abhyankar     }
8902b4ed282SShri Abhyankar     count++;
8912b4ed282SShri Abhyankar   }
8922b4ed282SShri Abhyankar   theend2:
8932b4ed282SShri Abhyankar   /* Optional user-defined check for line search step validity */
8940661309fSPeter Brune   if (snes->ops->postcheckstep) {
8950661309fSPeter Brune     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8962b4ed282SShri Abhyankar     if (changed_y) {
8972b4ed282SShri Abhyankar       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8987fe79bc4SShri Abhyankar       ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr);
8992b4ed282SShri Abhyankar     }
9002b4ed282SShri Abhyankar     if (changed_y || changed_w) { /* recompute the function if the step has changed */
9012b4ed282SShri Abhyankar       ierr = SNESComputeFunction(snes,w,g);
9022b4ed282SShri Abhyankar       if (snes->domainerror) {
9032b4ed282SShri Abhyankar         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
9042b4ed282SShri Abhyankar         PetscFunctionReturn(0);
9052b4ed282SShri Abhyankar       }
906*c2fc9fa9SBarry Smith       if (snes->ops->solve != SNESSolve_VISS) {
9077fe79bc4SShri Abhyankar         ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr);
9087fe79bc4SShri Abhyankar       } else {
9097fe79bc4SShri Abhyankar         ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
9107fe79bc4SShri Abhyankar       }
9117fe79bc4SShri Abhyankar 
9122b4ed282SShri Abhyankar       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
9132b4ed282SShri Abhyankar       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
91462cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
9152b4ed282SShri Abhyankar     }
9162b4ed282SShri Abhyankar   }
9172b4ed282SShri Abhyankar   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
9182b4ed282SShri Abhyankar   PetscFunctionReturn(0);
9192b4ed282SShri Abhyankar }
9202b4ed282SShri Abhyankar 
9212b4ed282SShri Abhyankar 
9225559b345SBarry Smith #undef __FUNCT__
9235559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds"
9245559b345SBarry Smith /*@
9252b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
9262b4ed282SShri Abhyankar 
9272b4ed282SShri Abhyankar    Input Parameters:
9282b4ed282SShri Abhyankar .  snes - the SNES context.
9292b4ed282SShri Abhyankar .  xl   - lower bound.
9302b4ed282SShri Abhyankar .  xu   - upper bound.
9312b4ed282SShri Abhyankar 
9322b4ed282SShri Abhyankar    Notes:
9332b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
93401f0b76bSBarry Smith    SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp().
93584c105d7SBarry Smith 
9362bd2b0e6SSatish Balay    Level: advanced
9372bd2b0e6SSatish Balay 
9385559b345SBarry Smith @*/
93969c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
9402b4ed282SShri Abhyankar {
941d923200aSJungho Lee   PetscErrorCode    ierr;
9426fd67740SBarry Smith   const PetscScalar *xxl,*xxu;
9436fd67740SBarry Smith   PetscInt          i,n, cnt = 0;
9442b4ed282SShri Abhyankar 
9452b4ed282SShri Abhyankar   PetscFunctionBegin;
9462b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
9472b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
9482b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
949a63bb30eSJed Brown   ierr = SNESGetFunction(snes,&snes->vec_func,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
950a63bb30eSJed Brown   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
9512b4ed282SShri 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);
9522b4ed282SShri 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);
953*c2fc9fa9SBarry Smith   ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);
9542176524cSBarry Smith   ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
9552176524cSBarry Smith   ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
956*c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xl);CHKERRQ(ierr);
957*c2fc9fa9SBarry Smith   ierr = VecDestroy(&snes->xu);CHKERRQ(ierr);
958*c2fc9fa9SBarry Smith   snes->xl = xl;
959*c2fc9fa9SBarry Smith   snes->xu = xu;
9606fd67740SBarry Smith   ierr = VecGetLocalSize(xl,&n);CHKERRQ(ierr);
9616fd67740SBarry Smith   ierr = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr);
9626fd67740SBarry Smith   ierr = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr);
9636fd67740SBarry Smith   for (i=0; i<n; i++) {
9646fd67740SBarry Smith     cnt += ((xxl[i] != SNES_VI_NINF) || (xxu[i] != SNES_VI_INF));
9656fd67740SBarry Smith   }
966*c2fc9fa9SBarry Smith   ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
9676fd67740SBarry Smith   ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr);
9686fd67740SBarry Smith   ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr);
9692b4ed282SShri Abhyankar   PetscFunctionReturn(0);
9702b4ed282SShri Abhyankar }
971693ddf92SShri Abhyankar 
97292c02d66SPeter Brune 
97392c02d66SPeter Brune EXTERN_C_BEGIN
97492c02d66SPeter Brune #undef __FUNCT__
97592c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_VI"
97692c02d66SPeter Brune PetscErrorCode  SNESLineSearchSetType_VI(SNES snes, SNESLineSearchType type)
97792c02d66SPeter Brune {
97892c02d66SPeter Brune   PetscErrorCode ierr;
97992c02d66SPeter Brune   PetscFunctionBegin;
98092c02d66SPeter Brune 
98192c02d66SPeter Brune   switch (type) {
98292c02d66SPeter Brune   case SNES_LS_BASIC:
98392c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr);
98492c02d66SPeter Brune     break;
98592c02d66SPeter Brune   case SNES_LS_BASIC_NONORMS:
98692c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr);
98792c02d66SPeter Brune     break;
98892c02d66SPeter Brune   case SNES_LS_QUADRATIC:
98992c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr);
99092c02d66SPeter Brune     break;
99192c02d66SPeter Brune   case SNES_LS_CUBIC:
99292c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr);
99392c02d66SPeter Brune     break;
99492c02d66SPeter Brune   default:
99592c02d66SPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
99692c02d66SPeter Brune     break;
99792c02d66SPeter Brune   }
99892c02d66SPeter Brune   snes->ls_type = type;
99992c02d66SPeter Brune   PetscFunctionReturn(0);
100092c02d66SPeter Brune }
100192c02d66SPeter Brune EXTERN_C_END
100292c02d66SPeter Brune 
10032b4ed282SShri Abhyankar #undef __FUNCT__
1004*c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VI"
1005*c2fc9fa9SBarry Smith PetscErrorCode SNESSetFromOptions_VI(SNES snes)
10062b4ed282SShri Abhyankar {
10072b4ed282SShri Abhyankar   PetscErrorCode  ierr;
1008*c2fc9fa9SBarry Smith   PetscBool       flg;
10092b4ed282SShri Abhyankar 
10102b4ed282SShri Abhyankar   PetscFunctionBegin;
1011*c2fc9fa9SBarry Smith   ierr = PetscOptionsHead("SNES VI options");CHKERRQ(ierr);
1012*c2fc9fa9SBarry Smith   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
1013*c2fc9fa9SBarry Smith   if (flg) {
1014*c2fc9fa9SBarry Smith     ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr);
10152b4ed282SShri Abhyankar   }
1016*c2fc9fa9SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1017ed70e96aSJungho Lee   PetscFunctionReturn(0);
1018ed70e96aSJungho Lee }
1019