xref: /petsc/src/snes/impls/vi/vi.c (revision b2566f29c2b6470df769aa9f7deb9e2726b0959e)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h>  /*I "petscsnes.h" I*/
21e25c274SJed Brown #include <petscdm.h>
3d0af7cd3SBarry Smith 
4d0af7cd3SBarry Smith #undef __FUNCT__
52176524cSBarry Smith #define __FUNCT__ "SNESVISetComputeVariableBounds"
62176524cSBarry Smith /*@C
72176524cSBarry Smith    SNESVISetComputeVariableBounds - Sets a function that is called to compute the variable bounds
82176524cSBarry Smith 
92176524cSBarry Smith    Input parameter
102176524cSBarry Smith +  snes - the SNES context
112176524cSBarry Smith -  compute - computes the bounds
122176524cSBarry Smith 
132bd2b0e6SSatish Balay    Level: advanced
142176524cSBarry Smith 
15c1c3a0ecSBarry Smith .seealso:   SNESVISetVariableBounds()
16c1c3a0ecSBarry Smith 
17aab9d709SJed Brown @*/
1877cdaf69SJed Brown PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec))
192176524cSBarry Smith {
2061589011SJed Brown   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec));
212176524cSBarry Smith 
222176524cSBarry Smith   PetscFunctionBegin;
2361589011SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
240005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",&f);CHKERRQ(ierr);
2563cdc2bdSPatrick Farrell   if (!f) {
2663cdc2bdSPatrick Farrell     ierr = SNESVISetComputeVariableBounds_VI(snes,compute);CHKERRQ(ierr);
2763cdc2bdSPatrick Farrell   } else {
2861589011SJed Brown     ierr = PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec)),(snes,compute));CHKERRQ(ierr);
2963cdc2bdSPatrick Farrell   }
302176524cSBarry Smith   PetscFunctionReturn(0);
312176524cSBarry Smith }
322176524cSBarry Smith 
3361589011SJed Brown #undef __FUNCT__
3461589011SJed Brown #define __FUNCT__ "SNESVISetComputeVariableBounds_VI"
3561589011SJed Brown PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute)
3661589011SJed Brown {
3761589011SJed Brown   PetscFunctionBegin;
3861589011SJed Brown   snes->ops->computevariablebounds = compute;
3961589011SJed Brown   PetscFunctionReturn(0);
4061589011SJed Brown }
412176524cSBarry Smith 
423c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/
432b4ed282SShri Abhyankar 
449308c137SBarry Smith #undef __FUNCT__
45ffdf2a20SBarry Smith #define __FUNCT__ "SNESVIMonitorResidual"
46ffdf2a20SBarry Smith PetscErrorCode  SNESVIMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
47ffdf2a20SBarry Smith {
48ffdf2a20SBarry Smith   PetscErrorCode ierr;
49ffdf2a20SBarry Smith   Vec            X, F, Finactive;
50ffdf2a20SBarry Smith   IS             isactive;
51ffdf2a20SBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
52ffdf2a20SBarry Smith 
53ffdf2a20SBarry Smith   PetscFunctionBegin;
54ffdf2a20SBarry Smith   ierr = SNESGetFunction(snes,&F,0,0);CHKERRQ(ierr);
55ffdf2a20SBarry Smith   ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr);
5687e98922SBarry Smith   ierr = SNESVIGetActiveSetIS(snes,X,F,&isactive);CHKERRQ(ierr);
57ffdf2a20SBarry Smith   ierr = VecDuplicate(F,&Finactive);CHKERRQ(ierr);
58ffdf2a20SBarry Smith   ierr = VecCopy(F,Finactive);CHKERRQ(ierr);
5987e98922SBarry Smith   ierr = VecISSet(Finactive,isactive,0.0);CHKERRQ(ierr);
60de34d3e9SBarry Smith   ierr = ISDestroy(&isactive);CHKERRQ(ierr);
6187e98922SBarry Smith   ierr = VecView(Finactive,viewer);CHKERRQ(ierr);
62ffdf2a20SBarry Smith   ierr = VecDestroy(&Finactive);CHKERRQ(ierr);
63ffdf2a20SBarry Smith   PetscFunctionReturn(0);
64ffdf2a20SBarry Smith }
65ffdf2a20SBarry Smith 
66ffdf2a20SBarry Smith #undef __FUNCT__
679308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI"
687087cfbeSBarry Smith PetscErrorCode  SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
699308c137SBarry Smith {
709308c137SBarry Smith   PetscErrorCode    ierr;
714d4332d5SBarry Smith   PetscViewer       viewer = (PetscViewer) dummy;
729308c137SBarry Smith   const PetscScalar *x,*xl,*xu,*f;
736fd67740SBarry Smith   PetscInt          i,n,act[2] = {0,0},fact[2],N;
746a9e2d46SJungho Lee   /* Number of components that actually hit the bounds (c.f. active variables) */
756a9e2d46SJungho Lee   PetscInt          act_bound[2] = {0,0},fact_bound[2];
76349187a7SBarry Smith   PetscReal         rnorm,fnorm,zerotolerance = snes->vizerotolerance;
779d1809e2SSatish Balay   double            tmp;
789308c137SBarry Smith 
799308c137SBarry Smith   PetscFunctionBegin;
804d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
819308c137SBarry Smith   ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr);
826fd67740SBarry Smith   ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr);
83c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
84c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
859308c137SBarry Smith   ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
863f731af5SBarry Smith   ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
879308c137SBarry Smith 
889308c137SBarry Smith   rnorm = 0.0;
899308c137SBarry Smith   for (i=0; i<n; i++) {
90349187a7SBarry Smith     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
91349187a7SBarry Smith     else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++;
92349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++;
93ce94432eSBarry Smith     else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Can never get here");
949308c137SBarry Smith   }
956a9e2d46SJungho Lee 
966a9e2d46SJungho Lee   for (i=0; i<n; i++) {
97349187a7SBarry Smith     if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++;
98349187a7SBarry Smith     else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++;
996a9e2d46SJungho Lee   }
1003f731af5SBarry Smith   ierr  = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr);
101c2fc9fa9SBarry Smith   ierr  = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
102c2fc9fa9SBarry Smith   ierr  = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
1039308c137SBarry Smith   ierr  = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr);
104*b2566f29SBarry Smith   ierr  = MPIU_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
105*b2566f29SBarry Smith   ierr  = MPIU_Allreduce(act,fact,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
106*b2566f29SBarry Smith   ierr  = MPIU_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
107f137e44dSBarry Smith   fnorm = PetscSqrtReal(fnorm);
1086fd67740SBarry Smith 
109649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1109d1809e2SSatish Balay   if (snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds);
1119d1809e2SSatish Balay   else tmp = 0.0;
11262151390SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES VI Function norm %g Active lower constraints %D/%D upper constraints %D/%D Percent of total %g Percent of bounded %g\n",its,(double)fnorm,fact[0],fact_bound[0],fact[1],fact_bound[1],((double)(fact[0]+fact[1]))/((double)N),tmp);CHKERRQ(ierr);
1136a9e2d46SJungho Lee 
114649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
1159308c137SBarry Smith   PetscFunctionReturn(0);
1169308c137SBarry Smith }
1179308c137SBarry Smith 
1182b4ed282SShri Abhyankar /*
1192b4ed282SShri Abhyankar      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1202b4ed282SShri 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
1212b4ed282SShri Abhyankar     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1222b4ed282SShri Abhyankar     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
1232b4ed282SShri Abhyankar */
1242b4ed282SShri Abhyankar #undef __FUNCT__
1252b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private"
126ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin)
1272b4ed282SShri Abhyankar {
1282b4ed282SShri Abhyankar   PetscReal      a1;
1292b4ed282SShri Abhyankar   PetscErrorCode ierr;
130ace3abfcSBarry Smith   PetscBool      hastranspose;
1312b4ed282SShri Abhyankar 
1322b4ed282SShri Abhyankar   PetscFunctionBegin;
1332b4ed282SShri Abhyankar   *ismin = PETSC_FALSE;
1342b4ed282SShri Abhyankar   ierr   = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1352b4ed282SShri Abhyankar   if (hastranspose) {
1362b4ed282SShri Abhyankar     /* Compute || J^T F|| */
1372b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
1382b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
1394839bfe8SBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
1402b4ed282SShri Abhyankar     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
1412b4ed282SShri Abhyankar   } else {
1422b4ed282SShri Abhyankar     Vec         work;
1432b4ed282SShri Abhyankar     PetscScalar result;
1442b4ed282SShri Abhyankar     PetscReal   wnorm;
1452b4ed282SShri Abhyankar 
1460298fd71SBarry Smith     ierr = VecSetRandom(W,NULL);CHKERRQ(ierr);
1472b4ed282SShri Abhyankar     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
1482b4ed282SShri Abhyankar     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
1492b4ed282SShri Abhyankar     ierr = MatMult(A,W,work);CHKERRQ(ierr);
1502b4ed282SShri Abhyankar     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
1516bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
1522b4ed282SShri Abhyankar     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
1534839bfe8SBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
1542b4ed282SShri Abhyankar     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
1552b4ed282SShri Abhyankar   }
1562b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1572b4ed282SShri Abhyankar }
1582b4ed282SShri Abhyankar 
1592b4ed282SShri Abhyankar /*
1602b4ed282SShri Abhyankar      Checks if J^T(F - J*X) = 0
1612b4ed282SShri Abhyankar */
1622b4ed282SShri Abhyankar #undef __FUNCT__
1632b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private"
1642b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
1652b4ed282SShri Abhyankar {
1662b4ed282SShri Abhyankar   PetscReal      a1,a2;
1672b4ed282SShri Abhyankar   PetscErrorCode ierr;
168ace3abfcSBarry Smith   PetscBool      hastranspose;
1692b4ed282SShri Abhyankar 
1702b4ed282SShri Abhyankar   PetscFunctionBegin;
1712b4ed282SShri Abhyankar   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
1722b4ed282SShri Abhyankar   if (hastranspose) {
1732b4ed282SShri Abhyankar     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
1742b4ed282SShri Abhyankar     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
1752b4ed282SShri Abhyankar 
1762b4ed282SShri Abhyankar     /* Compute || J^T W|| */
1772b4ed282SShri Abhyankar     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
1782b4ed282SShri Abhyankar     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
1792b4ed282SShri Abhyankar     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
1802b4ed282SShri Abhyankar     if (a1 != 0.0) {
1814839bfe8SBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
1822b4ed282SShri Abhyankar     }
1832b4ed282SShri Abhyankar   }
1842b4ed282SShri Abhyankar   PetscFunctionReturn(0);
1852b4ed282SShri Abhyankar }
1862b4ed282SShri Abhyankar 
1872b4ed282SShri Abhyankar /*
1888d359177SBarry Smith   SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm.
1892b4ed282SShri Abhyankar 
1902b4ed282SShri Abhyankar   Notes:
1912b4ed282SShri Abhyankar   The convergence criterion currently implemented is
1922b4ed282SShri Abhyankar   merit < abstol
1932b4ed282SShri Abhyankar   merit < rtol*merit_initial
1942b4ed282SShri Abhyankar */
1952b4ed282SShri Abhyankar #undef __FUNCT__
1968d359177SBarry Smith #define __FUNCT__ "SNESConvergedDefault_VI"
1978d359177SBarry Smith PetscErrorCode SNESConvergedDefault_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
1982b4ed282SShri Abhyankar {
1992b4ed282SShri Abhyankar   PetscErrorCode ierr;
2002b4ed282SShri Abhyankar 
2012b4ed282SShri Abhyankar   PetscFunctionBegin;
2022b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2032b4ed282SShri Abhyankar   PetscValidPointer(reason,6);
2042b4ed282SShri Abhyankar 
2052b4ed282SShri Abhyankar   *reason = SNES_CONVERGED_ITERATING;
2062b4ed282SShri Abhyankar 
2072b4ed282SShri Abhyankar   if (!it) {
2082b4ed282SShri Abhyankar     /* set parameter for default relative tolerance convergence test */
2097fe79bc4SShri Abhyankar     snes->ttol = fnorm*snes->rtol;
2102b4ed282SShri Abhyankar   }
2117fe79bc4SShri Abhyankar   if (fnorm != fnorm) {
2122b4ed282SShri Abhyankar     ierr    = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
2132b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FNORM_NAN;
2147fe79bc4SShri Abhyankar   } else if (fnorm < snes->abstol) {
2154839bfe8SBarry Smith     ierr    = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
2162b4ed282SShri Abhyankar     *reason = SNES_CONVERGED_FNORM_ABS;
2172b4ed282SShri Abhyankar   } else if (snes->nfuncs >= snes->max_funcs) {
2182b4ed282SShri Abhyankar     ierr    = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
2192b4ed282SShri Abhyankar     *reason = SNES_DIVERGED_FUNCTION_COUNT;
2202b4ed282SShri Abhyankar   }
2212b4ed282SShri Abhyankar 
2222b4ed282SShri Abhyankar   if (it && !*reason) {
2237fe79bc4SShri Abhyankar     if (fnorm < snes->ttol) {
2244839bfe8SBarry Smith       ierr    = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
2252b4ed282SShri Abhyankar       *reason = SNES_CONVERGED_FNORM_RELATIVE;
2262b4ed282SShri Abhyankar     }
2272b4ed282SShri Abhyankar   }
2282b4ed282SShri Abhyankar   PetscFunctionReturn(0);
2292b4ed282SShri Abhyankar }
2302b4ed282SShri Abhyankar 
231ee657d29SShri Abhyankar 
232c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */
233c1a5e521SShri Abhyankar /*
234c1a5e521SShri Abhyankar    SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n.
235c1a5e521SShri Abhyankar 
236c1a5e521SShri Abhyankar    Input Parameters:
237c1a5e521SShri Abhyankar .  SNES - nonlinear solver context
238c1a5e521SShri Abhyankar 
239c1a5e521SShri Abhyankar    Output Parameters:
240c1a5e521SShri Abhyankar .  X - Bound projected X
241c1a5e521SShri Abhyankar 
242c1a5e521SShri Abhyankar */
243c1a5e521SShri Abhyankar 
244c1a5e521SShri Abhyankar #undef __FUNCT__
245c1a5e521SShri Abhyankar #define __FUNCT__ "SNESVIProjectOntoBounds"
246c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X)
247c1a5e521SShri Abhyankar {
248c1a5e521SShri Abhyankar   PetscErrorCode    ierr;
249c1a5e521SShri Abhyankar   const PetscScalar *xl,*xu;
250c1a5e521SShri Abhyankar   PetscScalar       *x;
251c1a5e521SShri Abhyankar   PetscInt          i,n;
252c1a5e521SShri Abhyankar 
253c1a5e521SShri Abhyankar   PetscFunctionBegin;
254c1a5e521SShri Abhyankar   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
255c1a5e521SShri Abhyankar   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
256c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
257c2fc9fa9SBarry Smith   ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
258c1a5e521SShri Abhyankar 
259c1a5e521SShri Abhyankar   for (i = 0; i<n; i++) {
26010f5ae6bSBarry Smith     if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i];
26110f5ae6bSBarry Smith     else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i];
262c1a5e521SShri Abhyankar   }
263c1a5e521SShri Abhyankar   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
264c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
265c2fc9fa9SBarry Smith   ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
266c1a5e521SShri Abhyankar   PetscFunctionReturn(0);
267c1a5e521SShri Abhyankar }
268c1a5e521SShri Abhyankar 
26969c03620SShri Abhyankar 
27069c03620SShri Abhyankar #undef __FUNCT__
271693ddf92SShri Abhyankar #define __FUNCT__ "SNESVIGetActiveSetIS"
272693ddf92SShri Abhyankar /*
273693ddf92SShri Abhyankar    SNESVIGetActiveSetIndices - Gets the global indices for the active set variables
274693ddf92SShri Abhyankar 
275693ddf92SShri Abhyankar    Input parameter
276693ddf92SShri Abhyankar .  snes - the SNES context
277693ddf92SShri Abhyankar .  X    - the snes solution vector
278693ddf92SShri Abhyankar .  F    - the nonlinear function vector
279693ddf92SShri Abhyankar 
280693ddf92SShri Abhyankar    Output parameter
281693ddf92SShri Abhyankar .  ISact - active set index set
282693ddf92SShri Abhyankar  */
283693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS *ISact)
284d950fb63SShri Abhyankar {
285d950fb63SShri Abhyankar   PetscErrorCode    ierr;
286c2fc9fa9SBarry Smith   Vec               Xl=snes->xl,Xu=snes->xu;
287693ddf92SShri Abhyankar   const PetscScalar *x,*f,*xl,*xu;
288693ddf92SShri Abhyankar   PetscInt          *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0;
289349187a7SBarry Smith   PetscReal         zerotolerance = snes->vizerotolerance;
290d950fb63SShri Abhyankar 
291d950fb63SShri Abhyankar   PetscFunctionBegin;
292d950fb63SShri Abhyankar   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
293d950fb63SShri Abhyankar   ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr);
294fe835674SShri Abhyankar   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
295fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr);
296fe835674SShri Abhyankar   ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr);
297fe835674SShri Abhyankar   ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr);
298693ddf92SShri Abhyankar   /* Compute active set size */
299d950fb63SShri Abhyankar   for (i=0; i < nlocal;i++) {
300349187a7SBarry Smith     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) nloc_isact++;
301d950fb63SShri Abhyankar   }
302693ddf92SShri Abhyankar 
303785e854fSJed Brown   ierr = PetscMalloc1(nloc_isact,&idx_act);CHKERRQ(ierr);
304d950fb63SShri Abhyankar 
305693ddf92SShri Abhyankar   /* Set active set indices */
306d950fb63SShri Abhyankar   for (i=0; i < nlocal; i++) {
307349187a7SBarry Smith     if (!((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) idx_act[i1++] = ilow+i;
308d950fb63SShri Abhyankar   }
309d950fb63SShri Abhyankar 
310693ddf92SShri Abhyankar   /* Create active set IS */
311ce94432eSBarry Smith   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes),nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr);
312d950fb63SShri Abhyankar 
313fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
314fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr);
315fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr);
316fe835674SShri Abhyankar   ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
317d950fb63SShri Abhyankar   PetscFunctionReturn(0);
318d950fb63SShri Abhyankar }
319d950fb63SShri Abhyankar 
320693ddf92SShri Abhyankar #undef __FUNCT__
321693ddf92SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS"
322693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact)
323693ddf92SShri Abhyankar {
324693ddf92SShri Abhyankar   PetscErrorCode ierr;
325077aedafSJed Brown   PetscInt       rstart,rend;
326693ddf92SShri Abhyankar 
327693ddf92SShri Abhyankar   PetscFunctionBegin;
328693ddf92SShri Abhyankar   ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr);
329077aedafSJed Brown   ierr = VecGetOwnershipRange(X,&rstart,&rend);CHKERRQ(ierr);
330077aedafSJed Brown   ierr = ISComplement(*ISact,rstart,rend,ISinact);CHKERRQ(ierr);
331693ddf92SShri Abhyankar   PetscFunctionReturn(0);
332693ddf92SShri Abhyankar }
333693ddf92SShri Abhyankar 
334fe835674SShri Abhyankar #undef __FUNCT__
335fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm"
33610f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm)
337fe835674SShri Abhyankar {
338fe835674SShri Abhyankar   PetscErrorCode    ierr;
339fe835674SShri Abhyankar   const PetscScalar *x,*xl,*xu,*f;
340fe835674SShri Abhyankar   PetscInt          i,n;
341349187a7SBarry Smith   PetscReal         rnorm,zerotolerance = snes->vizerotolerance;;
342fe835674SShri Abhyankar 
343fe835674SShri Abhyankar   PetscFunctionBegin;
344fe835674SShri Abhyankar   ierr  = VecGetLocalSize(X,&n);CHKERRQ(ierr);
345c2fc9fa9SBarry Smith   ierr  = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr);
346c2fc9fa9SBarry Smith   ierr  = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr);
347fe835674SShri Abhyankar   ierr  = VecGetArrayRead(X,&x);CHKERRQ(ierr);
348fe835674SShri Abhyankar   ierr  = VecGetArrayRead(F,&f);CHKERRQ(ierr);
349fe835674SShri Abhyankar   rnorm = 0.0;
350fe835674SShri Abhyankar   for (i=0; i<n; i++) {
351349187a7SBarry Smith     if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + zerotolerance || (PetscRealPart(f[i]) <= 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - zerotolerance) || PetscRealPart(f[i]) >= 0.0))) rnorm += PetscRealPart(PetscConj(f[i])*f[i]);
3528f918527SKarl Rupp   }
353fe835674SShri Abhyankar   ierr   = VecRestoreArrayRead(F,&f);CHKERRQ(ierr);
354c2fc9fa9SBarry Smith   ierr   = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr);
355c2fc9fa9SBarry Smith   ierr   = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr);
356fe835674SShri Abhyankar   ierr   = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
357*b2566f29SBarry Smith   ierr   = MPIU_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
35862d1f40fSBarry Smith   *fnorm = PetscSqrtReal(*fnorm);
359fe835674SShri Abhyankar   PetscFunctionReturn(0);
360fe835674SShri Abhyankar }
361fe835674SShri Abhyankar 
36208da532bSDmitry Karpeev #undef __FUNCT__
36308da532bSDmitry Karpeev #define __FUNCT__ "SNESVIDMComputeVariableBounds"
36408da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu)
36508da532bSDmitry Karpeev {
36608da532bSDmitry Karpeev   PetscErrorCode ierr;
3676e111a19SKarl Rupp 
36808da532bSDmitry Karpeev   PetscFunctionBegin;
36908da532bSDmitry Karpeev   ierr = DMComputeVariableBounds(snes->dm, xl, xu);CHKERRQ(ierr);
37008da532bSDmitry Karpeev   PetscFunctionReturn(0);
37108da532bSDmitry Karpeev }
37208da532bSDmitry Karpeev 
3732f63a38cSShri Abhyankar 
3742b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
3752b4ed282SShri Abhyankar /*
376c2fc9fa9SBarry Smith    SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up
3772b4ed282SShri Abhyankar    of the SNESVI nonlinear solver.
3782b4ed282SShri Abhyankar 
3792b4ed282SShri Abhyankar    Input Parameter:
3802b4ed282SShri Abhyankar .  snes - the SNES context
3812b4ed282SShri Abhyankar 
3822b4ed282SShri Abhyankar    Application Interface Routine: SNESSetUp()
3832b4ed282SShri Abhyankar 
3842b4ed282SShri Abhyankar    Notes:
3852b4ed282SShri Abhyankar    For basic use of the SNES solvers, the user need not explicitly call
3862b4ed282SShri Abhyankar    SNESSetUp(), since these actions will automatically occur during
3872b4ed282SShri Abhyankar    the call to SNESSolve().
3882b4ed282SShri Abhyankar  */
3892b4ed282SShri Abhyankar #undef __FUNCT__
3902b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI"
3912b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes)
3922b4ed282SShri Abhyankar {
3932b4ed282SShri Abhyankar   PetscErrorCode ierr;
3942b4ed282SShri Abhyankar   PetscInt       i_start[3],i_end[3];
3952b4ed282SShri Abhyankar 
3962b4ed282SShri Abhyankar   PetscFunctionBegin;
397fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,1);CHKERRQ(ierr);
3986cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
3992b4ed282SShri Abhyankar 
40008da532bSDmitry Karpeev   if (!snes->ops->computevariablebounds && snes->dm) {
401a201590fSDmitry Karpeev     PetscBool flag;
402a201590fSDmitry Karpeev     ierr = DMHasVariableBounds(snes->dm, &flag);CHKERRQ(ierr);
4031aa26658SKarl Rupp 
40408da532bSDmitry Karpeev     snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds;
40508da532bSDmitry Karpeev   }
406a201590fSDmitry Karpeev   if (!snes->usersetbounds) {
407c2fc9fa9SBarry Smith     if (snes->ops->computevariablebounds) {
408c2fc9fa9SBarry Smith       if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);}
409c2fc9fa9SBarry Smith       if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);}
410c2fc9fa9SBarry Smith       ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr);
4111aa26658SKarl Rupp     } else if (!snes->xl && !snes->xu) {
4122176524cSBarry Smith       /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */
413c2fc9fa9SBarry Smith       ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr);
414e270355aSBarry Smith       ierr = VecSet(snes->xl,PETSC_NINFINITY);CHKERRQ(ierr);
415c2fc9fa9SBarry Smith       ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr);
416e270355aSBarry Smith       ierr = VecSet(snes->xu,PETSC_INFINITY);CHKERRQ(ierr);
4172b4ed282SShri Abhyankar     } else {
4182b4ed282SShri Abhyankar       /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */
4192b4ed282SShri Abhyankar       ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr);
420c2fc9fa9SBarry Smith       ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr);
421c2fc9fa9SBarry Smith       ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr);
4222b4ed282SShri 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]))
4232b4ed282SShri 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.");
4242b4ed282SShri Abhyankar     }
425a201590fSDmitry Karpeev   }
4262b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4272b4ed282SShri Abhyankar }
4282b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */
4292176524cSBarry Smith #undef __FUNCT__
4302176524cSBarry Smith #define __FUNCT__ "SNESReset_VI"
4312176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes)
4322176524cSBarry Smith {
4332176524cSBarry Smith   PetscErrorCode ierr;
4342176524cSBarry Smith 
4352176524cSBarry Smith   PetscFunctionBegin;
436c2fc9fa9SBarry Smith   ierr                = VecDestroy(&snes->xl);CHKERRQ(ierr);
437c2fc9fa9SBarry Smith   ierr                = VecDestroy(&snes->xu);CHKERRQ(ierr);
4382d6615e8SDmitry Karpeev   snes->usersetbounds = PETSC_FALSE;
4392176524cSBarry Smith   PetscFunctionReturn(0);
4402176524cSBarry Smith }
4412176524cSBarry Smith 
4422b4ed282SShri Abhyankar /*
4432b4ed282SShri Abhyankar    SNESDestroy_VI - Destroys the private SNES_VI context that was created
4442b4ed282SShri Abhyankar    with SNESCreate_VI().
4452b4ed282SShri Abhyankar 
4462b4ed282SShri Abhyankar    Input Parameter:
4472b4ed282SShri Abhyankar .  snes - the SNES context
4482b4ed282SShri Abhyankar 
4492b4ed282SShri Abhyankar    Application Interface Routine: SNESDestroy()
4502b4ed282SShri Abhyankar  */
4512b4ed282SShri Abhyankar #undef __FUNCT__
4522b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI"
4532b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes)
4542b4ed282SShri Abhyankar {
4552b4ed282SShri Abhyankar   PetscErrorCode ierr;
4562b4ed282SShri Abhyankar 
4572b4ed282SShri Abhyankar   PetscFunctionBegin;
4582b4ed282SShri Abhyankar   ierr = PetscFree(snes->data);CHKERRQ(ierr);
4592b4ed282SShri Abhyankar 
4602b4ed282SShri Abhyankar   /* clear composed functions */
461bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSet_C",NULL);CHKERRQ(ierr);
462bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSetMonitor_C",NULL);CHKERRQ(ierr);
4632b4ed282SShri Abhyankar   PetscFunctionReturn(0);
4642b4ed282SShri Abhyankar }
4657fe79bc4SShri Abhyankar 
4665559b345SBarry Smith #undef __FUNCT__
4675559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds"
4685559b345SBarry Smith /*@
4692b4ed282SShri Abhyankar    SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu.
4702b4ed282SShri Abhyankar 
4712b4ed282SShri Abhyankar    Input Parameters:
4722b4ed282SShri Abhyankar .  snes - the SNES context.
4732b4ed282SShri Abhyankar .  xl   - lower bound.
4742b4ed282SShri Abhyankar .  xu   - upper bound.
4752b4ed282SShri Abhyankar 
4762b4ed282SShri Abhyankar    Notes:
4772b4ed282SShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
47829eed3a4SBarry Smith    PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp().
47984c105d7SBarry Smith 
4802bd2b0e6SSatish Balay    Level: advanced
4812bd2b0e6SSatish Balay 
4825559b345SBarry Smith @*/
48369c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu)
4842b4ed282SShri Abhyankar {
48561589011SJed Brown   PetscErrorCode ierr,(*f)(SNES,Vec,Vec);
4862b4ed282SShri Abhyankar 
4872b4ed282SShri Abhyankar   PetscFunctionBegin;
4882b4ed282SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4892b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
4902b4ed282SShri Abhyankar   PetscValidHeaderSpecific(xu,VEC_CLASSID,3);
4910005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",&f);CHKERRQ(ierr);
49256e5e3feSPatrick Farrell   if (!f) {
49356e5e3feSPatrick Farrell     ierr = SNESVISetVariableBounds_VI(snes, xl, xu);CHKERRQ(ierr);
49456e5e3feSPatrick Farrell   } else {
49561589011SJed Brown     ierr = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr);
49656e5e3feSPatrick Farrell   }
497a201590fSDmitry Karpeev   snes->usersetbounds = PETSC_TRUE;
49861589011SJed Brown   PetscFunctionReturn(0);
49961589011SJed Brown }
50061589011SJed Brown 
50161589011SJed Brown #undef __FUNCT__
50261589011SJed Brown #define __FUNCT__ "SNESVISetVariableBounds_VI"
50361589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu)
50461589011SJed Brown {
50561589011SJed Brown   PetscErrorCode    ierr;
50661589011SJed Brown   const PetscScalar *xxl,*xxu;
50761589011SJed Brown   PetscInt          i,n, cnt = 0;
50861589011SJed Brown 
50961589011SJed Brown   PetscFunctionBegin;
5100298fd71SBarry Smith   ierr = SNESGetFunction(snes,&snes->vec_func,NULL,NULL);CHKERRQ(ierr);
511a63bb30eSJed Brown   if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
512077aedafSJed Brown   {
513077aedafSJed Brown     PetscInt xlN,xuN,N;
514077aedafSJed Brown     ierr = VecGetSize(xl,&xlN);CHKERRQ(ierr);
515077aedafSJed Brown     ierr = VecGetSize(xu,&xuN);CHKERRQ(ierr);
516077aedafSJed Brown     ierr = VecGetSize(snes->vec_func,&N);CHKERRQ(ierr);
517077aedafSJed Brown     if (xlN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xlN,N);
518077aedafSJed Brown     if (xuN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xuN,N);
519077aedafSJed Brown   }
5202176524cSBarry Smith   ierr     = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr);
5212176524cSBarry Smith   ierr     = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr);
522c2fc9fa9SBarry Smith   ierr     = VecDestroy(&snes->xl);CHKERRQ(ierr);
523c2fc9fa9SBarry Smith   ierr     = VecDestroy(&snes->xu);CHKERRQ(ierr);
524c2fc9fa9SBarry Smith   snes->xl = xl;
525c2fc9fa9SBarry Smith   snes->xu = xu;
5266fd67740SBarry Smith   ierr     = VecGetLocalSize(xl,&n);CHKERRQ(ierr);
5276fd67740SBarry Smith   ierr     = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr);
5286fd67740SBarry Smith   ierr     = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr);
529e270355aSBarry Smith   for (i=0; i<n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY));
5301aa26658SKarl Rupp 
531*b2566f29SBarry Smith   ierr = MPIU_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
5326fd67740SBarry Smith   ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr);
5336fd67740SBarry Smith   ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr);
5342b4ed282SShri Abhyankar   PetscFunctionReturn(0);
5352b4ed282SShri Abhyankar }
53692c02d66SPeter Brune 
5372b4ed282SShri Abhyankar #undef __FUNCT__
538c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VI"
5398c34d3f5SBarry Smith PetscErrorCode SNESSetFromOptions_VI(PetscOptions *PetscOptionsObject,SNES snes)
5402b4ed282SShri Abhyankar {
5412b4ed282SShri Abhyankar   PetscErrorCode ierr;
5428afaa268SBarry Smith   PetscBool      flg = PETSC_FALSE;
543639ea3faSPeter Brune   SNESLineSearch linesearch;
5442b4ed282SShri Abhyankar 
5452b4ed282SShri Abhyankar   PetscFunctionBegin;
546e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES VI options");CHKERRQ(ierr);
547349187a7SBarry Smith   ierr = PetscOptionsReal("-snes_vi_zero_tolerance","Tolerance for considering x[] value to be on a bound","None",snes->vizerotolerance,&snes->vizerotolerance,NULL);CHKERRQ(ierr);
548ffdf2a20SBarry Smith   ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","SNESMonitorResidual",flg,&flg,NULL);CHKERRQ(ierr);
549c2fc9fa9SBarry Smith   if (flg) {
5504d4332d5SBarry Smith     ierr = SNESMonitorSet(snes,SNESMonitorVI,PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),NULL);CHKERRQ(ierr);
5512b4ed282SShri Abhyankar   }
552de34d3e9SBarry Smith   flg = PETSC_FALSE;
553ffdf2a20SBarry Smith   ierr = PetscOptionsBool("-snes_vi_monitor_residual","Monitor residual all non-active variables; using zero for active constraints","SNESMonitorVIResidual",flg,&flg,NULL);CHKERRQ(ierr);
554ffdf2a20SBarry Smith   if (flg) {
555ffdf2a20SBarry Smith     ierr = SNESMonitorSet(snes,SNESVIMonitorResidual,PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)),NULL);CHKERRQ(ierr);
556ffdf2a20SBarry Smith   }
557639ea3faSPeter Brune   if (!snes->linesearch) {
5587601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
559639ea3faSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr);
560639ea3faSPeter Brune     ierr = SNESLineSearchBTSetAlpha(linesearch, 0.0);CHKERRQ(ierr);
561639ea3faSPeter Brune   }
562c2fc9fa9SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
563ed70e96aSJungho Lee   PetscFunctionReturn(0);
564ed70e96aSJungho Lee }
565