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); 25*63cdc2bdSPatrick Farrell if (!f) { 26*63cdc2bdSPatrick Farrell ierr = SNESVISetComputeVariableBounds_VI(snes,compute);CHKERRQ(ierr); 27*63cdc2bdSPatrick Farrell } else { 2861589011SJed Brown ierr = PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec)),(snes,compute));CHKERRQ(ierr); 29*63cdc2bdSPatrick 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 422176524cSBarry Smith #undef __FUNCT__ 43ed70e96aSJungho Lee #define __FUNCT__ "SNESVIComputeInactiveSetIS" 44d0af7cd3SBarry Smith /* 45ed70e96aSJungho Lee SNESVIComputeInactiveSetIS - Gets the global indices for the bogus inactive set variables 46d0af7cd3SBarry Smith 47d0af7cd3SBarry Smith Input parameter 48d0af7cd3SBarry Smith . snes - the SNES context 49d0af7cd3SBarry Smith . X - the snes solution vector 50d0af7cd3SBarry Smith 51d0af7cd3SBarry Smith Output parameter 52d0af7cd3SBarry Smith . ISact - active set index set 53d0af7cd3SBarry Smith 54d0af7cd3SBarry Smith */ 55ed70e96aSJungho Lee PetscErrorCode SNESVIComputeInactiveSetIS(Vec upper,Vec lower,Vec X,Vec F,IS *inact) 56d0af7cd3SBarry Smith { 57d0af7cd3SBarry Smith PetscErrorCode ierr; 58d0af7cd3SBarry Smith const PetscScalar *x,*xl,*xu,*f; 59d0af7cd3SBarry Smith PetscInt *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0; 60d0af7cd3SBarry Smith 61d0af7cd3SBarry Smith PetscFunctionBegin; 62d0af7cd3SBarry Smith ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 63d0af7cd3SBarry Smith ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr); 64d0af7cd3SBarry Smith ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 65d0af7cd3SBarry Smith ierr = VecGetArrayRead(lower,&xl);CHKERRQ(ierr); 66d0af7cd3SBarry Smith ierr = VecGetArrayRead(upper,&xu);CHKERRQ(ierr); 67d0af7cd3SBarry Smith ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 68d0af7cd3SBarry Smith /* Compute inactive set size */ 69d0af7cd3SBarry Smith for (i=0; i < nlocal; i++) { 70d0af7cd3SBarry Smith if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) nloc_isact++; 71d0af7cd3SBarry Smith } 72d0af7cd3SBarry Smith 73785e854fSJed Brown ierr = PetscMalloc1(nloc_isact,&idx_act);CHKERRQ(ierr); 74d0af7cd3SBarry Smith 75d0af7cd3SBarry Smith /* Set inactive set indices */ 76d0af7cd3SBarry Smith for (i=0; i < nlocal; i++) { 77d0af7cd3SBarry Smith if (((PetscRealPart(x[i]) > PetscRealPart(xl[i]) + 1.e-8 || (PetscRealPart(f[i]) < 0.0)) && ((PetscRealPart(x[i]) < PetscRealPart(xu[i]) - 1.e-8) || PetscRealPart(f[i]) > 0.0))) idx_act[i1++] = ilow+i; 78d0af7cd3SBarry Smith } 79d0af7cd3SBarry Smith 80d0af7cd3SBarry Smith /* Create inactive set IS */ 81ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)upper),nloc_isact,idx_act,PETSC_OWN_POINTER,inact);CHKERRQ(ierr); 82d0af7cd3SBarry Smith 83d0af7cd3SBarry Smith ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 84d0af7cd3SBarry Smith ierr = VecRestoreArrayRead(lower,&xl);CHKERRQ(ierr); 85d0af7cd3SBarry Smith ierr = VecRestoreArrayRead(upper,&xu);CHKERRQ(ierr); 86d0af7cd3SBarry Smith ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 87d0af7cd3SBarry Smith PetscFunctionReturn(0); 88d0af7cd3SBarry Smith } 89d0af7cd3SBarry Smith 903c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/ 912b4ed282SShri Abhyankar 929308c137SBarry Smith #undef __FUNCT__ 93ffdf2a20SBarry Smith #define __FUNCT__ "SNESVIMonitorResidual" 94ffdf2a20SBarry Smith PetscErrorCode SNESVIMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 95ffdf2a20SBarry Smith { 96ffdf2a20SBarry Smith PetscErrorCode ierr; 97ffdf2a20SBarry Smith Vec X, F, Finactive; 98ffdf2a20SBarry Smith IS isactive; 99ffdf2a20SBarry Smith PetscViewer viewer = (PetscViewer) dummy; 100ffdf2a20SBarry Smith 101ffdf2a20SBarry Smith PetscFunctionBegin; 102ffdf2a20SBarry Smith ierr = SNESGetFunction(snes,&F,0,0);CHKERRQ(ierr); 103ffdf2a20SBarry Smith ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr); 10487e98922SBarry Smith ierr = SNESVIGetActiveSetIS(snes,X,F,&isactive);CHKERRQ(ierr); 105ffdf2a20SBarry Smith ierr = VecDuplicate(F,&Finactive);CHKERRQ(ierr); 106ffdf2a20SBarry Smith ierr = VecCopy(F,Finactive);CHKERRQ(ierr); 10787e98922SBarry Smith ierr = VecISSet(Finactive,isactive,0.0);CHKERRQ(ierr); 108de34d3e9SBarry Smith ierr = ISDestroy(&isactive);CHKERRQ(ierr); 109ffdf2a20SBarry Smith if (!viewer) { 110ffdf2a20SBarry Smith viewer = PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)); 111ffdf2a20SBarry Smith } 11287e98922SBarry Smith ierr = VecView(Finactive,viewer);CHKERRQ(ierr); 113ffdf2a20SBarry Smith ierr = VecDestroy(&Finactive);CHKERRQ(ierr); 114ffdf2a20SBarry Smith PetscFunctionReturn(0); 115ffdf2a20SBarry Smith } 116ffdf2a20SBarry Smith 117ffdf2a20SBarry Smith #undef __FUNCT__ 1189308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI" 1197087cfbeSBarry Smith PetscErrorCode SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 1209308c137SBarry Smith { 1219308c137SBarry Smith PetscErrorCode ierr; 122ce94432eSBarry Smith PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); 1239308c137SBarry Smith const PetscScalar *x,*xl,*xu,*f; 1246fd67740SBarry Smith PetscInt i,n,act[2] = {0,0},fact[2],N; 1256a9e2d46SJungho Lee /* Number of components that actually hit the bounds (c.f. active variables) */ 1266a9e2d46SJungho Lee PetscInt act_bound[2] = {0,0},fact_bound[2]; 1279308c137SBarry Smith PetscReal rnorm,fnorm; 1289d1809e2SSatish Balay double tmp; 1299308c137SBarry Smith 1309308c137SBarry Smith PetscFunctionBegin; 1319308c137SBarry Smith ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr); 1326fd67740SBarry Smith ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr); 133c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 134c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 1359308c137SBarry Smith ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 1363f731af5SBarry Smith ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr); 1379308c137SBarry Smith 1389308c137SBarry Smith rnorm = 0.0; 1399308c137SBarry Smith for (i=0; i<n; i++) { 14010f5ae6bSBarry 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]); 141e400df7aSBarry Smith else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++; 142e400df7aSBarry Smith else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++; 143ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Can never get here"); 1449308c137SBarry Smith } 1456a9e2d46SJungho Lee 1466a9e2d46SJungho Lee for (i=0; i<n; i++) { 1476a9e2d46SJungho Lee if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++; 1486a9e2d46SJungho Lee else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++; 1496a9e2d46SJungho Lee } 1503f731af5SBarry Smith ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr); 151c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 152c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 1539308c137SBarry Smith ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 154ce94432eSBarry Smith ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 155ce94432eSBarry Smith ierr = MPI_Allreduce(act,fact,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 156ce94432eSBarry Smith ierr = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 157f137e44dSBarry Smith fnorm = PetscSqrtReal(fnorm); 1586fd67740SBarry Smith 159649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 1609d1809e2SSatish Balay if (snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds); 1619d1809e2SSatish Balay else tmp = 0.0; 1629d1809e2SSatish 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); 1636a9e2d46SJungho Lee 164649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 1659308c137SBarry Smith PetscFunctionReturn(0); 1669308c137SBarry Smith } 1679308c137SBarry Smith 1682b4ed282SShri Abhyankar /* 1692b4ed282SShri Abhyankar Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 1702b4ed282SShri 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 1712b4ed282SShri Abhyankar 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 1722b4ed282SShri Abhyankar for this trick. One assumes that the probability that W is in the null space of J is very, very small. 1732b4ed282SShri Abhyankar */ 1742b4ed282SShri Abhyankar #undef __FUNCT__ 1752b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private" 176ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 1772b4ed282SShri Abhyankar { 1782b4ed282SShri Abhyankar PetscReal a1; 1792b4ed282SShri Abhyankar PetscErrorCode ierr; 180ace3abfcSBarry Smith PetscBool hastranspose; 1812b4ed282SShri Abhyankar 1822b4ed282SShri Abhyankar PetscFunctionBegin; 1832b4ed282SShri Abhyankar *ismin = PETSC_FALSE; 1842b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 1852b4ed282SShri Abhyankar if (hastranspose) { 1862b4ed282SShri Abhyankar /* Compute || J^T F|| */ 1872b4ed282SShri Abhyankar ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 1882b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 1894839bfe8SBarry Smith ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr); 1902b4ed282SShri Abhyankar if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 1912b4ed282SShri Abhyankar } else { 1922b4ed282SShri Abhyankar Vec work; 1932b4ed282SShri Abhyankar PetscScalar result; 1942b4ed282SShri Abhyankar PetscReal wnorm; 1952b4ed282SShri Abhyankar 1960298fd71SBarry Smith ierr = VecSetRandom(W,NULL);CHKERRQ(ierr); 1972b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 1982b4ed282SShri Abhyankar ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 1992b4ed282SShri Abhyankar ierr = MatMult(A,W,work);CHKERRQ(ierr); 2002b4ed282SShri Abhyankar ierr = VecDot(F,work,&result);CHKERRQ(ierr); 2016bf464f9SBarry Smith ierr = VecDestroy(&work);CHKERRQ(ierr); 2022b4ed282SShri Abhyankar a1 = PetscAbsScalar(result)/(fnorm*wnorm); 2034839bfe8SBarry Smith ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr); 2042b4ed282SShri Abhyankar if (a1 < 1.e-4) *ismin = PETSC_TRUE; 2052b4ed282SShri Abhyankar } 2062b4ed282SShri Abhyankar PetscFunctionReturn(0); 2072b4ed282SShri Abhyankar } 2082b4ed282SShri Abhyankar 2092b4ed282SShri Abhyankar /* 2102b4ed282SShri Abhyankar Checks if J^T(F - J*X) = 0 2112b4ed282SShri Abhyankar */ 2122b4ed282SShri Abhyankar #undef __FUNCT__ 2132b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private" 2142b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 2152b4ed282SShri Abhyankar { 2162b4ed282SShri Abhyankar PetscReal a1,a2; 2172b4ed282SShri Abhyankar PetscErrorCode ierr; 218ace3abfcSBarry Smith PetscBool hastranspose; 2192b4ed282SShri Abhyankar 2202b4ed282SShri Abhyankar PetscFunctionBegin; 2212b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 2222b4ed282SShri Abhyankar if (hastranspose) { 2232b4ed282SShri Abhyankar ierr = MatMult(A,X,W1);CHKERRQ(ierr); 2242b4ed282SShri Abhyankar ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 2252b4ed282SShri Abhyankar 2262b4ed282SShri Abhyankar /* Compute || J^T W|| */ 2272b4ed282SShri Abhyankar ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 2282b4ed282SShri Abhyankar ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 2292b4ed282SShri Abhyankar ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 2302b4ed282SShri Abhyankar if (a1 != 0.0) { 2314839bfe8SBarry Smith ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr); 2322b4ed282SShri Abhyankar } 2332b4ed282SShri Abhyankar } 2342b4ed282SShri Abhyankar PetscFunctionReturn(0); 2352b4ed282SShri Abhyankar } 2362b4ed282SShri Abhyankar 2372b4ed282SShri Abhyankar /* 2388d359177SBarry Smith SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm. 2392b4ed282SShri Abhyankar 2402b4ed282SShri Abhyankar Notes: 2412b4ed282SShri Abhyankar The convergence criterion currently implemented is 2422b4ed282SShri Abhyankar merit < abstol 2432b4ed282SShri Abhyankar merit < rtol*merit_initial 2442b4ed282SShri Abhyankar */ 2452b4ed282SShri Abhyankar #undef __FUNCT__ 2468d359177SBarry Smith #define __FUNCT__ "SNESConvergedDefault_VI" 2478d359177SBarry Smith PetscErrorCode SNESConvergedDefault_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 2482b4ed282SShri Abhyankar { 2492b4ed282SShri Abhyankar PetscErrorCode ierr; 2502b4ed282SShri Abhyankar 2512b4ed282SShri Abhyankar PetscFunctionBegin; 2522b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2532b4ed282SShri Abhyankar PetscValidPointer(reason,6); 2542b4ed282SShri Abhyankar 2552b4ed282SShri Abhyankar *reason = SNES_CONVERGED_ITERATING; 2562b4ed282SShri Abhyankar 2572b4ed282SShri Abhyankar if (!it) { 2582b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 2597fe79bc4SShri Abhyankar snes->ttol = fnorm*snes->rtol; 2602b4ed282SShri Abhyankar } 2617fe79bc4SShri Abhyankar if (fnorm != fnorm) { 2622b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 2632b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FNORM_NAN; 2647fe79bc4SShri Abhyankar } else if (fnorm < snes->abstol) { 2654839bfe8SBarry Smith ierr = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr); 2662b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_ABS; 2672b4ed282SShri Abhyankar } else if (snes->nfuncs >= snes->max_funcs) { 2682b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 2692b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FUNCTION_COUNT; 2702b4ed282SShri Abhyankar } 2712b4ed282SShri Abhyankar 2722b4ed282SShri Abhyankar if (it && !*reason) { 2737fe79bc4SShri Abhyankar if (fnorm < snes->ttol) { 2744839bfe8SBarry Smith ierr = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr); 2752b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_RELATIVE; 2762b4ed282SShri Abhyankar } 2772b4ed282SShri Abhyankar } 2782b4ed282SShri Abhyankar PetscFunctionReturn(0); 2792b4ed282SShri Abhyankar } 2802b4ed282SShri Abhyankar 281ee657d29SShri Abhyankar 282c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */ 283c1a5e521SShri Abhyankar /* 284c1a5e521SShri Abhyankar SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n. 285c1a5e521SShri Abhyankar 286c1a5e521SShri Abhyankar Input Parameters: 287c1a5e521SShri Abhyankar . SNES - nonlinear solver context 288c1a5e521SShri Abhyankar 289c1a5e521SShri Abhyankar Output Parameters: 290c1a5e521SShri Abhyankar . X - Bound projected X 291c1a5e521SShri Abhyankar 292c1a5e521SShri Abhyankar */ 293c1a5e521SShri Abhyankar 294c1a5e521SShri Abhyankar #undef __FUNCT__ 295c1a5e521SShri Abhyankar #define __FUNCT__ "SNESVIProjectOntoBounds" 296c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X) 297c1a5e521SShri Abhyankar { 298c1a5e521SShri Abhyankar PetscErrorCode ierr; 299c1a5e521SShri Abhyankar const PetscScalar *xl,*xu; 300c1a5e521SShri Abhyankar PetscScalar *x; 301c1a5e521SShri Abhyankar PetscInt i,n; 302c1a5e521SShri Abhyankar 303c1a5e521SShri Abhyankar PetscFunctionBegin; 304c1a5e521SShri Abhyankar ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 305c1a5e521SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 306c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 307c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 308c1a5e521SShri Abhyankar 309c1a5e521SShri Abhyankar for (i = 0; i<n; i++) { 31010f5ae6bSBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i]; 31110f5ae6bSBarry Smith else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i]; 312c1a5e521SShri Abhyankar } 313c1a5e521SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 314c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 315c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 316c1a5e521SShri Abhyankar PetscFunctionReturn(0); 317c1a5e521SShri Abhyankar } 318c1a5e521SShri Abhyankar 31969c03620SShri Abhyankar 32069c03620SShri Abhyankar #undef __FUNCT__ 321693ddf92SShri Abhyankar #define __FUNCT__ "SNESVIGetActiveSetIS" 322693ddf92SShri Abhyankar /* 323693ddf92SShri Abhyankar SNESVIGetActiveSetIndices - Gets the global indices for the active set variables 324693ddf92SShri Abhyankar 325693ddf92SShri Abhyankar Input parameter 326693ddf92SShri Abhyankar . snes - the SNES context 327693ddf92SShri Abhyankar . X - the snes solution vector 328693ddf92SShri Abhyankar . F - the nonlinear function vector 329693ddf92SShri Abhyankar 330693ddf92SShri Abhyankar Output parameter 331693ddf92SShri Abhyankar . ISact - active set index set 332693ddf92SShri Abhyankar */ 333693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS *ISact) 334d950fb63SShri Abhyankar { 335d950fb63SShri Abhyankar PetscErrorCode ierr; 336c2fc9fa9SBarry Smith Vec Xl=snes->xl,Xu=snes->xu; 337693ddf92SShri Abhyankar const PetscScalar *x,*f,*xl,*xu; 338693ddf92SShri Abhyankar PetscInt *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0; 339d950fb63SShri Abhyankar 340d950fb63SShri Abhyankar PetscFunctionBegin; 341d950fb63SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 342d950fb63SShri Abhyankar ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr); 343fe835674SShri Abhyankar ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 344fe835674SShri Abhyankar ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr); 345fe835674SShri Abhyankar ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr); 346fe835674SShri Abhyankar ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 347693ddf92SShri Abhyankar /* Compute active set size */ 348d950fb63SShri Abhyankar for (i=0; i < nlocal;i++) { 34910f5ae6bSBarry 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++; 350d950fb63SShri Abhyankar } 351693ddf92SShri Abhyankar 352785e854fSJed Brown ierr = PetscMalloc1(nloc_isact,&idx_act);CHKERRQ(ierr); 353d950fb63SShri Abhyankar 354693ddf92SShri Abhyankar /* Set active set indices */ 355d950fb63SShri Abhyankar for (i=0; i < nlocal; i++) { 35610f5ae6bSBarry 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; 357d950fb63SShri Abhyankar } 358d950fb63SShri Abhyankar 359693ddf92SShri Abhyankar /* Create active set IS */ 360ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes),nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr); 361d950fb63SShri Abhyankar 362fe835674SShri Abhyankar ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 363fe835674SShri Abhyankar ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr); 364fe835674SShri Abhyankar ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr); 365fe835674SShri Abhyankar ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 366d950fb63SShri Abhyankar PetscFunctionReturn(0); 367d950fb63SShri Abhyankar } 368d950fb63SShri Abhyankar 369693ddf92SShri Abhyankar #undef __FUNCT__ 370693ddf92SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS" 371693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact) 372693ddf92SShri Abhyankar { 373693ddf92SShri Abhyankar PetscErrorCode ierr; 374077aedafSJed Brown PetscInt rstart,rend; 375693ddf92SShri Abhyankar 376693ddf92SShri Abhyankar PetscFunctionBegin; 377693ddf92SShri Abhyankar ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr); 378077aedafSJed Brown ierr = VecGetOwnershipRange(X,&rstart,&rend);CHKERRQ(ierr); 379077aedafSJed Brown ierr = ISComplement(*ISact,rstart,rend,ISinact);CHKERRQ(ierr); 380693ddf92SShri Abhyankar PetscFunctionReturn(0); 381693ddf92SShri Abhyankar } 382693ddf92SShri Abhyankar 383fe835674SShri Abhyankar #undef __FUNCT__ 384fe835674SShri Abhyankar #define __FUNCT__ "SNESVIComputeInactiveSetFnorm" 38510f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm) 386fe835674SShri Abhyankar { 387fe835674SShri Abhyankar PetscErrorCode ierr; 388fe835674SShri Abhyankar const PetscScalar *x,*xl,*xu,*f; 389fe835674SShri Abhyankar PetscInt i,n; 390fe835674SShri Abhyankar PetscReal rnorm; 391fe835674SShri Abhyankar 392fe835674SShri Abhyankar PetscFunctionBegin; 393fe835674SShri Abhyankar ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 394c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 395c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 396fe835674SShri Abhyankar ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 397fe835674SShri Abhyankar ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 398fe835674SShri Abhyankar rnorm = 0.0; 399fe835674SShri Abhyankar for (i=0; i<n; i++) { 40010f5ae6bSBarry 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]); 4018f918527SKarl Rupp } 402fe835674SShri Abhyankar ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 403c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 404c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 405fe835674SShri Abhyankar ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 406ce94432eSBarry Smith ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 40762d1f40fSBarry Smith *fnorm = PetscSqrtReal(*fnorm); 408fe835674SShri Abhyankar PetscFunctionReturn(0); 409fe835674SShri Abhyankar } 410fe835674SShri Abhyankar 41108da532bSDmitry Karpeev #undef __FUNCT__ 41208da532bSDmitry Karpeev #define __FUNCT__ "SNESVIDMComputeVariableBounds" 41308da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu) 41408da532bSDmitry Karpeev { 41508da532bSDmitry Karpeev PetscErrorCode ierr; 4166e111a19SKarl Rupp 41708da532bSDmitry Karpeev PetscFunctionBegin; 41808da532bSDmitry Karpeev ierr = DMComputeVariableBounds(snes->dm, xl, xu);CHKERRQ(ierr); 41908da532bSDmitry Karpeev PetscFunctionReturn(0); 42008da532bSDmitry Karpeev } 42108da532bSDmitry Karpeev 4222f63a38cSShri Abhyankar 4232b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 4242b4ed282SShri Abhyankar /* 425c2fc9fa9SBarry Smith SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up 4262b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 4272b4ed282SShri Abhyankar 4282b4ed282SShri Abhyankar Input Parameter: 4292b4ed282SShri Abhyankar . snes - the SNES context 4302b4ed282SShri Abhyankar 4312b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 4322b4ed282SShri Abhyankar 4332b4ed282SShri Abhyankar Notes: 4342b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 4352b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 4362b4ed282SShri Abhyankar the call to SNESSolve(). 4372b4ed282SShri Abhyankar */ 4382b4ed282SShri Abhyankar #undef __FUNCT__ 4392b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI" 4402b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes) 4412b4ed282SShri Abhyankar { 4422b4ed282SShri Abhyankar PetscErrorCode ierr; 4432b4ed282SShri Abhyankar PetscInt i_start[3],i_end[3]; 4442b4ed282SShri Abhyankar 4452b4ed282SShri Abhyankar PetscFunctionBegin; 446fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,1);CHKERRQ(ierr); 4476cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 4482b4ed282SShri Abhyankar 44908da532bSDmitry Karpeev if (!snes->ops->computevariablebounds && snes->dm) { 450a201590fSDmitry Karpeev PetscBool flag; 451a201590fSDmitry Karpeev ierr = DMHasVariableBounds(snes->dm, &flag);CHKERRQ(ierr); 4521aa26658SKarl Rupp 45308da532bSDmitry Karpeev snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds; 45408da532bSDmitry Karpeev } 455a201590fSDmitry Karpeev if (!snes->usersetbounds) { 456c2fc9fa9SBarry Smith if (snes->ops->computevariablebounds) { 457c2fc9fa9SBarry Smith if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);} 458c2fc9fa9SBarry Smith if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);} 459c2fc9fa9SBarry Smith ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr); 4601aa26658SKarl Rupp } else if (!snes->xl && !snes->xu) { 4612176524cSBarry Smith /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */ 462c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr); 463e270355aSBarry Smith ierr = VecSet(snes->xl,PETSC_NINFINITY);CHKERRQ(ierr); 464c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr); 465e270355aSBarry Smith ierr = VecSet(snes->xu,PETSC_INFINITY);CHKERRQ(ierr); 4662b4ed282SShri Abhyankar } else { 4672b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 4682b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 469c2fc9fa9SBarry Smith ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr); 470c2fc9fa9SBarry Smith ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr); 4712b4ed282SShri 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])) 4722b4ed282SShri 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."); 4732b4ed282SShri Abhyankar } 474a201590fSDmitry Karpeev } 4752b4ed282SShri Abhyankar PetscFunctionReturn(0); 4762b4ed282SShri Abhyankar } 4772b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 4782176524cSBarry Smith #undef __FUNCT__ 4792176524cSBarry Smith #define __FUNCT__ "SNESReset_VI" 4802176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes) 4812176524cSBarry Smith { 4822176524cSBarry Smith PetscErrorCode ierr; 4832176524cSBarry Smith 4842176524cSBarry Smith PetscFunctionBegin; 485c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xl);CHKERRQ(ierr); 486c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xu);CHKERRQ(ierr); 4872d6615e8SDmitry Karpeev snes->usersetbounds = PETSC_FALSE; 4882176524cSBarry Smith PetscFunctionReturn(0); 4892176524cSBarry Smith } 4902176524cSBarry Smith 4912b4ed282SShri Abhyankar /* 4922b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 4932b4ed282SShri Abhyankar with SNESCreate_VI(). 4942b4ed282SShri Abhyankar 4952b4ed282SShri Abhyankar Input Parameter: 4962b4ed282SShri Abhyankar . snes - the SNES context 4972b4ed282SShri Abhyankar 4982b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 4992b4ed282SShri Abhyankar */ 5002b4ed282SShri Abhyankar #undef __FUNCT__ 5012b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI" 5022b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes) 5032b4ed282SShri Abhyankar { 5042b4ed282SShri Abhyankar PetscErrorCode ierr; 5052b4ed282SShri Abhyankar 5062b4ed282SShri Abhyankar PetscFunctionBegin; 5072b4ed282SShri Abhyankar ierr = PetscFree(snes->data);CHKERRQ(ierr); 5082b4ed282SShri Abhyankar 5092b4ed282SShri Abhyankar /* clear composed functions */ 510bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSet_C",NULL);CHKERRQ(ierr); 511bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSetMonitor_C",NULL);CHKERRQ(ierr); 5122b4ed282SShri Abhyankar PetscFunctionReturn(0); 5132b4ed282SShri Abhyankar } 5147fe79bc4SShri Abhyankar 5155559b345SBarry Smith #undef __FUNCT__ 5165559b345SBarry Smith #define __FUNCT__ "SNESVISetVariableBounds" 5175559b345SBarry Smith /*@ 5182b4ed282SShri Abhyankar SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 5192b4ed282SShri Abhyankar 5202b4ed282SShri Abhyankar Input Parameters: 5212b4ed282SShri Abhyankar . snes - the SNES context. 5222b4ed282SShri Abhyankar . xl - lower bound. 5232b4ed282SShri Abhyankar . xu - upper bound. 5242b4ed282SShri Abhyankar 5252b4ed282SShri Abhyankar Notes: 5262b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 527e270355aSBarry Smith PETSC_INFINITY and PETSC_NINFINITY respectively during SNESSetUp(). 52884c105d7SBarry Smith 5292bd2b0e6SSatish Balay Level: advanced 5302bd2b0e6SSatish Balay 5315559b345SBarry Smith @*/ 53269c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 5332b4ed282SShri Abhyankar { 53461589011SJed Brown PetscErrorCode ierr,(*f)(SNES,Vec,Vec); 5352b4ed282SShri Abhyankar 5362b4ed282SShri Abhyankar PetscFunctionBegin; 5372b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5382b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 5392b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 5400005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",&f);CHKERRQ(ierr); 54156e5e3feSPatrick Farrell if (!f) { 54256e5e3feSPatrick Farrell ierr = SNESVISetVariableBounds_VI(snes, xl, xu);CHKERRQ(ierr); 54356e5e3feSPatrick Farrell } else { 54461589011SJed Brown ierr = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr); 54556e5e3feSPatrick Farrell } 546a201590fSDmitry Karpeev snes->usersetbounds = PETSC_TRUE; 54761589011SJed Brown PetscFunctionReturn(0); 54861589011SJed Brown } 54961589011SJed Brown 55061589011SJed Brown #undef __FUNCT__ 55161589011SJed Brown #define __FUNCT__ "SNESVISetVariableBounds_VI" 55261589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu) 55361589011SJed Brown { 55461589011SJed Brown PetscErrorCode ierr; 55561589011SJed Brown const PetscScalar *xxl,*xxu; 55661589011SJed Brown PetscInt i,n, cnt = 0; 55761589011SJed Brown 55861589011SJed Brown PetscFunctionBegin; 5590298fd71SBarry Smith ierr = SNESGetFunction(snes,&snes->vec_func,NULL,NULL);CHKERRQ(ierr); 560a63bb30eSJed Brown if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first"); 561077aedafSJed Brown { 562077aedafSJed Brown PetscInt xlN,xuN,N; 563077aedafSJed Brown ierr = VecGetSize(xl,&xlN);CHKERRQ(ierr); 564077aedafSJed Brown ierr = VecGetSize(xu,&xuN);CHKERRQ(ierr); 565077aedafSJed Brown ierr = VecGetSize(snes->vec_func,&N);CHKERRQ(ierr); 566077aedafSJed Brown if (xlN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xlN,N); 567077aedafSJed Brown if (xuN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xuN,N); 568077aedafSJed Brown } 5692176524cSBarry Smith ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr); 5702176524cSBarry Smith ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr); 571c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xl);CHKERRQ(ierr); 572c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xu);CHKERRQ(ierr); 573c2fc9fa9SBarry Smith snes->xl = xl; 574c2fc9fa9SBarry Smith snes->xu = xu; 5756fd67740SBarry Smith ierr = VecGetLocalSize(xl,&n);CHKERRQ(ierr); 5766fd67740SBarry Smith ierr = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr); 5776fd67740SBarry Smith ierr = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr); 578e270355aSBarry Smith for (i=0; i<n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY)); 5791aa26658SKarl Rupp 580ce94432eSBarry Smith ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 5816fd67740SBarry Smith ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr); 5826fd67740SBarry Smith ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr); 5832b4ed282SShri Abhyankar PetscFunctionReturn(0); 5842b4ed282SShri Abhyankar } 58592c02d66SPeter Brune 5862b4ed282SShri Abhyankar #undef __FUNCT__ 587c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VI" 5888c34d3f5SBarry Smith PetscErrorCode SNESSetFromOptions_VI(PetscOptions *PetscOptionsObject,SNES snes) 5892b4ed282SShri Abhyankar { 5902b4ed282SShri Abhyankar PetscErrorCode ierr; 5918afaa268SBarry Smith PetscBool flg = PETSC_FALSE; 592639ea3faSPeter Brune SNESLineSearch linesearch; 5932b4ed282SShri Abhyankar 5942b4ed282SShri Abhyankar PetscFunctionBegin; 595e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES VI options");CHKERRQ(ierr); 596ffdf2a20SBarry Smith ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","SNESMonitorResidual",flg,&flg,NULL);CHKERRQ(ierr); 597c2fc9fa9SBarry Smith if (flg) { 598c2fc9fa9SBarry Smith ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr); 5992b4ed282SShri Abhyankar } 600de34d3e9SBarry Smith flg = PETSC_FALSE; 601ffdf2a20SBarry Smith ierr = PetscOptionsBool("-snes_vi_monitor_residual","Monitor residual all non-active variables; using zero for active constraints","SNESMonitorVIResidual",flg,&flg,NULL);CHKERRQ(ierr); 602ffdf2a20SBarry Smith if (flg) { 603ffdf2a20SBarry Smith ierr = SNESMonitorSet(snes,SNESVIMonitorResidual,PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)),NULL);CHKERRQ(ierr); 604ffdf2a20SBarry Smith } 605639ea3faSPeter Brune if (!snes->linesearch) { 6067601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 607639ea3faSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr); 608639ea3faSPeter Brune ierr = SNESLineSearchBTSetAlpha(linesearch, 0.0);CHKERRQ(ierr); 609639ea3faSPeter Brune } 610c2fc9fa9SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 611ed70e96aSJungho Lee PetscFunctionReturn(0); 612ed70e96aSJungho Lee } 613