1639ea3faSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2d0af7cd3SBarry Smith 3d0af7cd3SBarry Smith #undef __FUNCT__ 42176524cSBarry Smith #define __FUNCT__ "SNESVISetComputeVariableBounds" 52176524cSBarry Smith /*@C 62176524cSBarry Smith SNESVISetComputeVariableBounds - Sets a function that is called to compute the variable bounds 72176524cSBarry Smith 82176524cSBarry Smith Input parameter 92176524cSBarry Smith + snes - the SNES context 102176524cSBarry Smith - compute - computes the bounds 112176524cSBarry Smith 122bd2b0e6SSatish Balay Level: advanced 132176524cSBarry Smith 14c1c3a0ecSBarry Smith .seealso: SNESVISetVariableBounds() 15c1c3a0ecSBarry Smith 16aab9d709SJed Brown @*/ 1777cdaf69SJed Brown PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec)) 182176524cSBarry Smith { 1961589011SJed Brown PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec)); 202176524cSBarry Smith 212176524cSBarry Smith PetscFunctionBegin; 2261589011SJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2361589011SJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr); 24f450aa47SBarry Smith if (!f) {ierr = SNESSetType(snes,SNESVINEWTONRSLS);CHKERRQ(ierr);} 2561589011SJed Brown ierr = PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec)),(snes,compute));CHKERRQ(ierr); 262176524cSBarry Smith PetscFunctionReturn(0); 272176524cSBarry Smith } 282176524cSBarry Smith 2961589011SJed Brown EXTERN_C_BEGIN 3061589011SJed Brown #undef __FUNCT__ 3161589011SJed Brown #define __FUNCT__ "SNESVISetComputeVariableBounds_VI" 3261589011SJed Brown PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute) 3361589011SJed Brown { 3461589011SJed Brown PetscFunctionBegin; 3561589011SJed Brown snes->ops->computevariablebounds = compute; 3661589011SJed Brown PetscFunctionReturn(0); 3761589011SJed Brown } 3861589011SJed Brown EXTERN_C_END 392176524cSBarry Smith 402176524cSBarry Smith #undef __FUNCT__ 41ed70e96aSJungho Lee #define __FUNCT__ "SNESVIComputeInactiveSetIS" 42d0af7cd3SBarry Smith /* 43ed70e96aSJungho Lee SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables 44d0af7cd3SBarry Smith 45d0af7cd3SBarry Smith Input parameter 46d0af7cd3SBarry Smith . snes - the SNES context 47d0af7cd3SBarry Smith . X - the snes solution vector 48d0af7cd3SBarry Smith 49d0af7cd3SBarry Smith Output parameter 50d0af7cd3SBarry Smith . ISact - active set index set 51d0af7cd3SBarry Smith 52d0af7cd3SBarry Smith */ 53ed70e96aSJungho Lee PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS *inact) 54d0af7cd3SBarry Smith { 55d0af7cd3SBarry Smith PetscErrorCode ierr; 56d0af7cd3SBarry Smith const PetscScalar *x,*xl,*xu,*f; 57d0af7cd3SBarry Smith PetscInt *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0; 58d0af7cd3SBarry Smith 59d0af7cd3SBarry Smith PetscFunctionBegin; 60d0af7cd3SBarry Smith ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 61d0af7cd3SBarry Smith ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr); 62d0af7cd3SBarry Smith ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 63d0af7cd3SBarry Smith ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr); 64d0af7cd3SBarry Smith ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr); 65d0af7cd3SBarry Smith ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 66d0af7cd3SBarry Smith /* Compute inactive set size */ 67d0af7cd3SBarry Smith for (i=0; i < nlocal; i++) { 68d0af7cd3SBarry 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++; 69d0af7cd3SBarry Smith } 70d0af7cd3SBarry Smith 71d0af7cd3SBarry Smith ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 72d0af7cd3SBarry Smith 73d0af7cd3SBarry Smith /* Set inactive set indices */ 74d0af7cd3SBarry Smith for (i=0; i < nlocal; i++) { 75d0af7cd3SBarry 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; 76d0af7cd3SBarry Smith } 77d0af7cd3SBarry Smith 78d0af7cd3SBarry Smith /* Create inactive set IS */ 79*ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)upper),nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr); 80d0af7cd3SBarry Smith 81d0af7cd3SBarry Smith ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 82d0af7cd3SBarry Smith ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr); 83d0af7cd3SBarry Smith ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr); 84d0af7cd3SBarry Smith ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 85d0af7cd3SBarry Smith PetscFunctionReturn(0); 86d0af7cd3SBarry Smith } 87d0af7cd3SBarry Smith 883c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/ 892b4ed282SShri Abhyankar 909308c137SBarry Smith #undef __FUNCT__ 919308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI" 927087cfbeSBarry Smith PetscErrorCode SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 939308c137SBarry Smith { 949308c137SBarry Smith PetscErrorCode ierr; 95*ce94432eSBarry Smith PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); 969308c137SBarry Smith const PetscScalar *x,*xl,*xu,*f; 976fd67740SBarry Smith PetscInt i,n,act[2] = {0,0},fact[2],N; 986a9e2d46SJungho Lee /* Number of components that actually hit the bounds (c.f. active variables) */ 996a9e2d46SJungho Lee PetscInt act_bound[2] = {0,0},fact_bound[2]; 1009308c137SBarry Smith PetscReal rnorm,fnorm; 1019d1809e2SSatish Balay double tmp; 1029308c137SBarry Smith 1039308c137SBarry Smith PetscFunctionBegin; 1049308c137SBarry Smith ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr); 1056fd67740SBarry Smith ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr); 106c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 107c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 1089308c137SBarry Smith ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 1093f731af5SBarry Smith ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr); 1109308c137SBarry Smith 1119308c137SBarry Smith rnorm = 0.0; 1129308c137SBarry Smith for (i=0; i<n; i++) { 11310f5ae6bSBarry 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]); 114e400df7aSBarry Smith else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++; 115e400df7aSBarry Smith else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++; 116*ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Can never get here"); 1179308c137SBarry Smith } 1186a9e2d46SJungho Lee 1196a9e2d46SJungho Lee for (i=0; i<n; i++) { 1206a9e2d46SJungho Lee if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++; 1216a9e2d46SJungho Lee else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++; 1226a9e2d46SJungho Lee } 1233f731af5SBarry Smith ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr); 124c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 125c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 1269308c137SBarry Smith ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 127*ce94432eSBarry Smith ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 128*ce94432eSBarry Smith ierr = MPI_Allreduce(act,fact,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 129*ce94432eSBarry Smith ierr = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 130f137e44dSBarry Smith fnorm = PetscSqrtReal(fnorm); 1316fd67740SBarry Smith 132649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 1339d1809e2SSatish Balay if (snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds); 1349d1809e2SSatish Balay else tmp = 0.0; 1359d1809e2SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES VI Function norm %14.12e Active lower constraints %D/%D upper constraints %D/%D Percent of total %g Percent of bounded %g\n",its,(double)fnorm,fact[0],fact_bound[0],fact[1],fact_bound[1],((double)(fact[0]+fact[1]))/((double)N),tmp);CHKERRQ(ierr); 1366a9e2d46SJungho Lee 137649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 1389308c137SBarry Smith PetscFunctionReturn(0); 1399308c137SBarry Smith } 1409308c137SBarry Smith 1412b4ed282SShri Abhyankar /* 1422b4ed282SShri Abhyankar Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 1432b4ed282SShri 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 1442b4ed282SShri Abhyankar 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 1452b4ed282SShri Abhyankar for this trick. One assumes that the probability that W is in the null space of J is very, very small. 1462b4ed282SShri Abhyankar */ 1472b4ed282SShri Abhyankar #undef __FUNCT__ 1482b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private" 149ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 1502b4ed282SShri Abhyankar { 1512b4ed282SShri Abhyankar PetscReal a1; 1522b4ed282SShri Abhyankar PetscErrorCode ierr; 153ace3abfcSBarry Smith PetscBool hastranspose; 1542b4ed282SShri Abhyankar 1552b4ed282SShri Abhyankar PetscFunctionBegin; 1562b4ed282SShri Abhyankar *ismin = PETSC_FALSE; 1572b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 1582b4ed282SShri Abhyankar if (hastranspose) { 1592b4ed282SShri Abhyankar /* Compute || J^T F|| */ 1602b4ed282SShri Abhyankar ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 1612b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 1624839bfe8SBarry Smith ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr); 1632b4ed282SShri Abhyankar if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 1642b4ed282SShri Abhyankar } else { 1652b4ed282SShri Abhyankar Vec work; 1662b4ed282SShri Abhyankar PetscScalar result; 1672b4ed282SShri Abhyankar PetscReal wnorm; 1682b4ed282SShri Abhyankar 1690298fd71SBarry Smith ierr = VecSetRandom(W,NULL);CHKERRQ(ierr); 1702b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 1712b4ed282SShri Abhyankar ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 1722b4ed282SShri Abhyankar ierr = MatMult(A,W,work);CHKERRQ(ierr); 1732b4ed282SShri Abhyankar ierr = VecDot(F,work,&result);CHKERRQ(ierr); 1746bf464f9SBarry Smith ierr = VecDestroy(&work);CHKERRQ(ierr); 1752b4ed282SShri Abhyankar a1 = PetscAbsScalar(result)/(fnorm*wnorm); 1764839bfe8SBarry Smith ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr); 1772b4ed282SShri Abhyankar if (a1 < 1.e-4) *ismin = PETSC_TRUE; 1782b4ed282SShri Abhyankar } 1792b4ed282SShri Abhyankar PetscFunctionReturn(0); 1802b4ed282SShri Abhyankar } 1812b4ed282SShri Abhyankar 1822b4ed282SShri Abhyankar /* 1832b4ed282SShri Abhyankar Checks if J^T(F - J*X) = 0 1842b4ed282SShri Abhyankar */ 1852b4ed282SShri Abhyankar #undef __FUNCT__ 1862b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private" 1872b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 1882b4ed282SShri Abhyankar { 1892b4ed282SShri Abhyankar PetscReal a1,a2; 1902b4ed282SShri Abhyankar PetscErrorCode ierr; 191ace3abfcSBarry Smith PetscBool hastranspose; 1922b4ed282SShri Abhyankar 1932b4ed282SShri Abhyankar PetscFunctionBegin; 1942b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 1952b4ed282SShri Abhyankar if (hastranspose) { 1962b4ed282SShri Abhyankar ierr = MatMult(A,X,W1);CHKERRQ(ierr); 1972b4ed282SShri Abhyankar ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 1982b4ed282SShri Abhyankar 1992b4ed282SShri Abhyankar /* Compute || J^T W|| */ 2002b4ed282SShri Abhyankar ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 2012b4ed282SShri Abhyankar ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 2022b4ed282SShri Abhyankar ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 2032b4ed282SShri Abhyankar if (a1 != 0.0) { 2044839bfe8SBarry Smith ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr); 2052b4ed282SShri Abhyankar } 2062b4ed282SShri Abhyankar } 2072b4ed282SShri Abhyankar PetscFunctionReturn(0); 2082b4ed282SShri Abhyankar } 2092b4ed282SShri Abhyankar 2102b4ed282SShri Abhyankar /* 2112b4ed282SShri Abhyankar SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm. 2122b4ed282SShri Abhyankar 2132b4ed282SShri Abhyankar Notes: 2142b4ed282SShri Abhyankar The convergence criterion currently implemented is 2152b4ed282SShri Abhyankar merit < abstol 2162b4ed282SShri Abhyankar merit < rtol*merit_initial 2172b4ed282SShri Abhyankar */ 2182b4ed282SShri Abhyankar #undef __FUNCT__ 2192b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI" 2207fe79bc4SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 2212b4ed282SShri Abhyankar { 2222b4ed282SShri Abhyankar PetscErrorCode ierr; 2232b4ed282SShri Abhyankar 2242b4ed282SShri Abhyankar PetscFunctionBegin; 2252b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2262b4ed282SShri Abhyankar PetscValidPointer(reason,6); 2272b4ed282SShri Abhyankar 2282b4ed282SShri Abhyankar *reason = SNES_CONVERGED_ITERATING; 2292b4ed282SShri Abhyankar 2302b4ed282SShri Abhyankar if (!it) { 2312b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 2327fe79bc4SShri Abhyankar snes->ttol = fnorm*snes->rtol; 2332b4ed282SShri Abhyankar } 2347fe79bc4SShri Abhyankar if (fnorm != fnorm) { 2352b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 2362b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FNORM_NAN; 2377fe79bc4SShri Abhyankar } else if (fnorm < snes->abstol) { 2384839bfe8SBarry Smith ierr = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr); 2392b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_ABS; 2402b4ed282SShri Abhyankar } else if (snes->nfuncs >= snes->max_funcs) { 2412b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 2422b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FUNCTION_COUNT; 2432b4ed282SShri Abhyankar } 2442b4ed282SShri Abhyankar 2452b4ed282SShri Abhyankar if (it && !*reason) { 2467fe79bc4SShri Abhyankar if (fnorm < snes->ttol) { 2474839bfe8SBarry Smith ierr = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr); 2482b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_RELATIVE; 2492b4ed282SShri Abhyankar } 2502b4ed282SShri Abhyankar } 2512b4ed282SShri Abhyankar PetscFunctionReturn(0); 2522b4ed282SShri Abhyankar } 2532b4ed282SShri Abhyankar 254ee657d29SShri Abhyankar 255c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */ 256c1a5e521SShri Abhyankar /* 257c1a5e521SShri Abhyankar SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n. 258c1a5e521SShri Abhyankar 259c1a5e521SShri Abhyankar Input Parameters: 260c1a5e521SShri Abhyankar . SNES - nonlinear solver context 261c1a5e521SShri Abhyankar 262c1a5e521SShri Abhyankar Output Parameters: 263c1a5e521SShri Abhyankar . X - Bound projected X 264c1a5e521SShri Abhyankar 265c1a5e521SShri Abhyankar */ 266c1a5e521SShri Abhyankar 267c1a5e521SShri Abhyankar #undef __FUNCT__ 268c1a5e521SShri Abhyankar #define __FUNCT__ "SNESVIProjectOntoBounds" 269c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X) 270c1a5e521SShri Abhyankar { 271c1a5e521SShri Abhyankar PetscErrorCode ierr; 272c1a5e521SShri Abhyankar const PetscScalar *xl,*xu; 273c1a5e521SShri Abhyankar PetscScalar *x; 274c1a5e521SShri Abhyankar PetscInt i,n; 275c1a5e521SShri Abhyankar 276c1a5e521SShri Abhyankar PetscFunctionBegin; 277c1a5e521SShri Abhyankar ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 278c1a5e521SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 279c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 280c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 281c1a5e521SShri Abhyankar 282c1a5e521SShri Abhyankar for (i = 0; i<n; i++) { 28310f5ae6bSBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i]; 28410f5ae6bSBarry Smith else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i]; 285c1a5e521SShri Abhyankar } 286c1a5e521SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 287c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 288c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 289c1a5e521SShri Abhyankar PetscFunctionReturn(0); 290c1a5e521SShri Abhyankar } 291c1a5e521SShri Abhyankar 29269c03620SShri Abhyankar 29369c03620SShri Abhyankar #undef __FUNCT__ 294693ddf92SShri Abhyankar #define __FUNCT__ "SNESVIGetActiveSetIS" 295693ddf92SShri Abhyankar /* 296693ddf92SShri Abhyankar SNESVIGetActiveSetIndices - Gets the global indices for the active set variables 297693ddf92SShri Abhyankar 298693ddf92SShri Abhyankar Input parameter 299693ddf92SShri Abhyankar . snes - the SNES context 300693ddf92SShri Abhyankar . X - the snes solution vector 301693ddf92SShri Abhyankar . F - the nonlinear function vector 302693ddf92SShri Abhyankar 303693ddf92SShri Abhyankar Output parameter 304693ddf92SShri Abhyankar . ISact - active set index set 305693ddf92SShri Abhyankar */ 306693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS *ISact) 307d950fb63SShri Abhyankar { 308d950fb63SShri Abhyankar PetscErrorCode ierr; 309c2fc9fa9SBarry Smith Vec Xl=snes->xl,Xu=snes->xu; 310693ddf92SShri Abhyankar const PetscScalar *x,*f,*xl,*xu; 311693ddf92SShri Abhyankar PetscInt *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0; 312d950fb63SShri Abhyankar 313d950fb63SShri Abhyankar PetscFunctionBegin; 314d950fb63SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 315d950fb63SShri Abhyankar ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr); 316fe835674SShri Abhyankar ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 317fe835674SShri Abhyankar ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr); 318fe835674SShri Abhyankar ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr); 319fe835674SShri Abhyankar ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 320693ddf92SShri Abhyankar /* Compute active set size */ 321d950fb63SShri Abhyankar for (i=0; i < nlocal;i++) { 32210f5ae6bSBarry 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++; 323d950fb63SShri Abhyankar } 324693ddf92SShri Abhyankar 325d950fb63SShri Abhyankar ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 326d950fb63SShri Abhyankar 327693ddf92SShri Abhyankar /* Set active set indices */ 328d950fb63SShri Abhyankar for (i=0; i < nlocal; i++) { 32910f5ae6bSBarry 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; 330d950fb63SShri Abhyankar } 331d950fb63SShri Abhyankar 332693ddf92SShri Abhyankar /* Create active set IS */ 333*ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes),nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr); 334d950fb63SShri Abhyankar 335fe835674SShri Abhyankar ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 336fe835674SShri Abhyankar ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr); 337fe835674SShri Abhyankar ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr); 338fe835674SShri Abhyankar ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 339d950fb63SShri Abhyankar PetscFunctionReturn(0); 340d950fb63SShri Abhyankar } 341d950fb63SShri Abhyankar 342693ddf92SShri Abhyankar #undef __FUNCT__ 343693ddf92SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS" 344693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact) 345693ddf92SShri Abhyankar { 346693ddf92SShri Abhyankar PetscErrorCode ierr; 347077aedafSJed Brown PetscInt rstart,rend; 348693ddf92SShri Abhyankar 349693ddf92SShri Abhyankar PetscFunctionBegin; 350693ddf92SShri Abhyankar ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr); 351077aedafSJed Brown ierr = VecGetOwnershipRange(X,&rstart,&rend);CHKERRQ(ierr); 352077aedafSJed Brown ierr = ISComplement(*ISact,rstart,rend,ISinact);CHKERRQ(ierr); 353693ddf92SShri Abhyankar PetscFunctionReturn(0); 354693ddf92SShri Abhyankar } 355693ddf92SShri Abhyankar 356fe835674SShri Abhyankar #undef __FUNCT__ 357fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm" 35810f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm) 359fe835674SShri Abhyankar { 360fe835674SShri Abhyankar PetscErrorCode ierr; 361fe835674SShri Abhyankar const PetscScalar *x,*xl,*xu,*f; 362fe835674SShri Abhyankar PetscInt i,n; 363fe835674SShri Abhyankar PetscReal rnorm; 364fe835674SShri Abhyankar 365fe835674SShri Abhyankar PetscFunctionBegin; 366fe835674SShri Abhyankar ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 367c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 368c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 369fe835674SShri Abhyankar ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 370fe835674SShri Abhyankar ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 371fe835674SShri Abhyankar rnorm = 0.0; 372fe835674SShri Abhyankar for (i=0; i<n; i++) { 37310f5ae6bSBarry 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]); 3748f918527SKarl Rupp } 375fe835674SShri Abhyankar ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 376c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 377c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 378fe835674SShri Abhyankar ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 379*ce94432eSBarry Smith ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 38062d1f40fSBarry Smith *fnorm = PetscSqrtReal(*fnorm); 381fe835674SShri Abhyankar PetscFunctionReturn(0); 382fe835674SShri Abhyankar } 383fe835674SShri Abhyankar 38408da532bSDmitry Karpeev #undef __FUNCT__ 38508da532bSDmitry Karpeev #define __FUNCT__ "SNESVIDMComputeVariableBounds" 38608da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu) 38708da532bSDmitry Karpeev { 38808da532bSDmitry Karpeev PetscErrorCode ierr; 3896e111a19SKarl Rupp 39008da532bSDmitry Karpeev PetscFunctionBegin; 39108da532bSDmitry Karpeev ierr = DMComputeVariableBounds(snes->dm, xl, xu);CHKERRQ(ierr); 39208da532bSDmitry Karpeev PetscFunctionReturn(0); 39308da532bSDmitry Karpeev } 39408da532bSDmitry Karpeev 3952f63a38cSShri Abhyankar 3962b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 3972b4ed282SShri Abhyankar /* 398c2fc9fa9SBarry Smith SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up 3992b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 4002b4ed282SShri Abhyankar 4012b4ed282SShri Abhyankar Input Parameter: 4022b4ed282SShri Abhyankar . snes - the SNES context 4032b4ed282SShri Abhyankar . x - the solution vector 4042b4ed282SShri Abhyankar 4052b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 4062b4ed282SShri Abhyankar 4072b4ed282SShri Abhyankar Notes: 4082b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 4092b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 4102b4ed282SShri Abhyankar the call to SNESSolve(). 4112b4ed282SShri Abhyankar */ 4122b4ed282SShri Abhyankar #undef __FUNCT__ 4132b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI" 4142b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes) 4152b4ed282SShri Abhyankar { 4162b4ed282SShri Abhyankar PetscErrorCode ierr; 4172b4ed282SShri Abhyankar PetscInt i_start[3],i_end[3]; 4182b4ed282SShri Abhyankar 4192b4ed282SShri Abhyankar PetscFunctionBegin; 4209bd66eb0SPeter Brune ierr = SNESDefaultGetWork(snes,1);CHKERRQ(ierr); 4216cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 4222b4ed282SShri Abhyankar 42308da532bSDmitry Karpeev if (!snes->ops->computevariablebounds && snes->dm) { 424a201590fSDmitry Karpeev PetscBool flag; 425a201590fSDmitry Karpeev ierr = DMHasVariableBounds(snes->dm, &flag);CHKERRQ(ierr); 4261aa26658SKarl Rupp 42708da532bSDmitry Karpeev snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds; 42808da532bSDmitry Karpeev } 429a201590fSDmitry Karpeev if (!snes->usersetbounds) { 430c2fc9fa9SBarry Smith if (snes->ops->computevariablebounds) { 431c2fc9fa9SBarry Smith if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);} 432c2fc9fa9SBarry Smith if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);} 433c2fc9fa9SBarry Smith ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr); 4341aa26658SKarl Rupp } else if (!snes->xl && !snes->xu) { 4352176524cSBarry Smith /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */ 436c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr); 437c2fc9fa9SBarry Smith ierr = VecSet(snes->xl,SNES_VI_NINF);CHKERRQ(ierr); 438c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr); 439c2fc9fa9SBarry Smith ierr = VecSet(snes->xu,SNES_VI_INF);CHKERRQ(ierr); 4402b4ed282SShri Abhyankar } else { 4412b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 4422b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 443c2fc9fa9SBarry Smith ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr); 444c2fc9fa9SBarry Smith ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr); 4452b4ed282SShri 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])) 4462b4ed282SShri 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."); 4472b4ed282SShri Abhyankar } 448a201590fSDmitry Karpeev } 4492b4ed282SShri Abhyankar PetscFunctionReturn(0); 4502b4ed282SShri Abhyankar } 4512b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 4522176524cSBarry Smith #undef __FUNCT__ 4532176524cSBarry Smith #define __FUNCT__ "SNESReset_VI" 4542176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes) 4552176524cSBarry Smith { 4562176524cSBarry Smith PetscErrorCode ierr; 4572176524cSBarry Smith 4582176524cSBarry Smith PetscFunctionBegin; 459c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xl);CHKERRQ(ierr); 460c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xu);CHKERRQ(ierr); 4612d6615e8SDmitry Karpeev snes->usersetbounds = PETSC_FALSE; 4622176524cSBarry Smith PetscFunctionReturn(0); 4632176524cSBarry Smith } 4642176524cSBarry Smith 4652b4ed282SShri Abhyankar /* 4662b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 4672b4ed282SShri Abhyankar with SNESCreate_VI(). 4682b4ed282SShri Abhyankar 4692b4ed282SShri Abhyankar Input Parameter: 4702b4ed282SShri Abhyankar . snes - the SNES context 4712b4ed282SShri Abhyankar 4722b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 4732b4ed282SShri Abhyankar */ 4742b4ed282SShri Abhyankar #undef __FUNCT__ 4752b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI" 4762b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes) 4772b4ed282SShri Abhyankar { 4782b4ed282SShri Abhyankar PetscErrorCode ierr; 4792b4ed282SShri Abhyankar 4802b4ed282SShri Abhyankar PetscFunctionBegin; 4812b4ed282SShri Abhyankar ierr = PetscFree(snes->data);CHKERRQ(ierr); 4822b4ed282SShri Abhyankar 4832b4ed282SShri Abhyankar /* clear composed functions */ 4840298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",NULL);CHKERRQ(ierr); 4850298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",NULL);CHKERRQ(ierr); 4862b4ed282SShri Abhyankar PetscFunctionReturn(0); 4872b4ed282SShri Abhyankar } 4887fe79bc4SShri Abhyankar 4895559b345SBarry Smith #undef __FUNCT__ 4905559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds" 4915559b345SBarry Smith /*@ 4922b4ed282SShri Abhyankar SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 4932b4ed282SShri Abhyankar 4942b4ed282SShri Abhyankar Input Parameters: 4952b4ed282SShri Abhyankar . snes - the SNES context. 4962b4ed282SShri Abhyankar . xl - lower bound. 4972b4ed282SShri Abhyankar . xu - upper bound. 4982b4ed282SShri Abhyankar 4992b4ed282SShri Abhyankar Notes: 5002b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 50101f0b76bSBarry Smith SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp(). 50284c105d7SBarry Smith 5032bd2b0e6SSatish Balay Level: advanced 5042bd2b0e6SSatish Balay 5055559b345SBarry Smith @*/ 50669c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 5072b4ed282SShri Abhyankar { 50861589011SJed Brown PetscErrorCode ierr,(*f)(SNES,Vec,Vec); 5092b4ed282SShri Abhyankar 5102b4ed282SShri Abhyankar PetscFunctionBegin; 5112b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5122b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 5132b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 51461589011SJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr); 515f450aa47SBarry Smith if (!f) {ierr = SNESSetType(snes,SNESVINEWTONRSLS);CHKERRQ(ierr);} 51661589011SJed Brown ierr = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr); 517a201590fSDmitry Karpeev snes->usersetbounds = PETSC_TRUE; 51861589011SJed Brown PetscFunctionReturn(0); 51961589011SJed Brown } 52061589011SJed Brown 52161589011SJed Brown EXTERN_C_BEGIN 52261589011SJed Brown #undef __FUNCT__ 52361589011SJed Brown #define __FUNCT__ "SNESVISetVariableBounds_VI" 52461589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu) 52561589011SJed Brown { 52661589011SJed Brown PetscErrorCode ierr; 52761589011SJed Brown const PetscScalar *xxl,*xxu; 52861589011SJed Brown PetscInt i,n, cnt = 0; 52961589011SJed Brown 53061589011SJed Brown PetscFunctionBegin; 5310298fd71SBarry Smith ierr = SNESGetFunction(snes,&snes->vec_func,NULL,NULL);CHKERRQ(ierr); 532a63bb30eSJed Brown if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first"); 533077aedafSJed Brown { 534077aedafSJed Brown PetscInt xlN,xuN,N; 535077aedafSJed Brown ierr = VecGetSize(xl,&xlN);CHKERRQ(ierr); 536077aedafSJed Brown ierr = VecGetSize(xu,&xuN);CHKERRQ(ierr); 537077aedafSJed Brown ierr = VecGetSize(snes->vec_func,&N);CHKERRQ(ierr); 538077aedafSJed Brown if (xlN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xlN,N); 539077aedafSJed Brown if (xuN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xuN,N); 540077aedafSJed Brown } 541f450aa47SBarry Smith ierr = SNESSetType(snes,SNESVINEWTONRSLS);CHKERRQ(ierr); 5422176524cSBarry Smith ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr); 5432176524cSBarry Smith ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr); 544c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xl);CHKERRQ(ierr); 545c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xu);CHKERRQ(ierr); 546c2fc9fa9SBarry Smith snes->xl = xl; 547c2fc9fa9SBarry Smith snes->xu = xu; 5486fd67740SBarry Smith ierr = VecGetLocalSize(xl,&n);CHKERRQ(ierr); 5496fd67740SBarry Smith ierr = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr); 5506fd67740SBarry Smith ierr = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr); 5511aa26658SKarl Rupp for (i=0; i<n; i++) cnt += ((xxl[i] != SNES_VI_NINF) || (xxu[i] != SNES_VI_INF)); 5521aa26658SKarl Rupp 553*ce94432eSBarry Smith ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 5546fd67740SBarry Smith ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr); 5556fd67740SBarry Smith ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr); 5562b4ed282SShri Abhyankar PetscFunctionReturn(0); 5572b4ed282SShri Abhyankar } 55861589011SJed Brown EXTERN_C_END 55992c02d66SPeter Brune 5602b4ed282SShri Abhyankar #undef __FUNCT__ 561c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VI" 562c2fc9fa9SBarry Smith PetscErrorCode SNESSetFromOptions_VI(SNES snes) 5632b4ed282SShri Abhyankar { 5642b4ed282SShri Abhyankar PetscErrorCode ierr; 565c2fc9fa9SBarry Smith PetscBool flg; 566639ea3faSPeter Brune SNESLineSearch linesearch; 5672b4ed282SShri Abhyankar 5682b4ed282SShri Abhyankar PetscFunctionBegin; 569c2fc9fa9SBarry Smith ierr = PetscOptionsHead("SNES VI options");CHKERRQ(ierr); 570c2fc9fa9SBarry Smith ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 571c2fc9fa9SBarry Smith if (flg) { 572c2fc9fa9SBarry Smith ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr); 5732b4ed282SShri Abhyankar } 574639ea3faSPeter Brune if (!snes->linesearch) { 575639ea3faSPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 576639ea3faSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr); 577639ea3faSPeter Brune ierr = SNESLineSearchBTSetAlpha(linesearch, 0.0);CHKERRQ(ierr); 578639ea3faSPeter Brune } 579c2fc9fa9SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 580ed70e96aSJungho Lee PetscFunctionReturn(0); 581ed70e96aSJungho Lee } 582