1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 21e25c274SJed Brown #include <petscdm.h> 3d0af7cd3SBarry Smith 42176524cSBarry Smith /*@C 52176524cSBarry Smith SNESVISetComputeVariableBounds - Sets a function that is called to compute the variable bounds 62176524cSBarry Smith 72176524cSBarry Smith Input parameter 82176524cSBarry Smith + snes - the SNES context 92176524cSBarry Smith - compute - computes the bounds 102176524cSBarry Smith 112bd2b0e6SSatish Balay Level: advanced 122176524cSBarry Smith 13c1c3a0ecSBarry Smith .seealso: SNESVISetVariableBounds() 14c1c3a0ecSBarry Smith 15aab9d709SJed Brown @*/ 1677cdaf69SJed Brown PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec)) 172176524cSBarry Smith { 1861589011SJed Brown PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec)); 192176524cSBarry Smith 202176524cSBarry Smith PetscFunctionBegin; 2161589011SJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 220005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",&f);CHKERRQ(ierr); 2363cdc2bdSPatrick Farrell if (!f) { 2463cdc2bdSPatrick Farrell ierr = SNESVISetComputeVariableBounds_VI(snes,compute);CHKERRQ(ierr); 2563cdc2bdSPatrick Farrell } else { 2661589011SJed Brown ierr = PetscUseMethod(snes,"SNESVISetComputeVariableBounds_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec)),(snes,compute));CHKERRQ(ierr); 2763cdc2bdSPatrick Farrell } 282176524cSBarry Smith PetscFunctionReturn(0); 292176524cSBarry Smith } 302176524cSBarry Smith 3161589011SJed Brown PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes,SNESVIComputeVariableBoundsFunction compute) 3261589011SJed Brown { 3361589011SJed Brown PetscFunctionBegin; 3461589011SJed Brown snes->ops->computevariablebounds = compute; 3561589011SJed Brown PetscFunctionReturn(0); 3661589011SJed Brown } 372176524cSBarry Smith 383c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/ 392b4ed282SShri Abhyankar 40ffdf2a20SBarry Smith PetscErrorCode SNESVIMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 41ffdf2a20SBarry Smith { 42ffdf2a20SBarry Smith PetscErrorCode ierr; 43ffdf2a20SBarry Smith Vec X, F, Finactive; 44ffdf2a20SBarry Smith IS isactive; 45ffdf2a20SBarry Smith PetscViewer viewer = (PetscViewer) dummy; 46ffdf2a20SBarry Smith 47ffdf2a20SBarry Smith PetscFunctionBegin; 48ffdf2a20SBarry Smith ierr = SNESGetFunction(snes,&F,0,0);CHKERRQ(ierr); 49ffdf2a20SBarry Smith ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr); 5087e98922SBarry Smith ierr = SNESVIGetActiveSetIS(snes,X,F,&isactive);CHKERRQ(ierr); 51ffdf2a20SBarry Smith ierr = VecDuplicate(F,&Finactive);CHKERRQ(ierr); 52ffdf2a20SBarry Smith ierr = VecCopy(F,Finactive);CHKERRQ(ierr); 5387e98922SBarry Smith ierr = VecISSet(Finactive,isactive,0.0);CHKERRQ(ierr); 54de34d3e9SBarry Smith ierr = ISDestroy(&isactive);CHKERRQ(ierr); 5587e98922SBarry Smith ierr = VecView(Finactive,viewer);CHKERRQ(ierr); 56ffdf2a20SBarry Smith ierr = VecDestroy(&Finactive);CHKERRQ(ierr); 57ffdf2a20SBarry Smith PetscFunctionReturn(0); 58ffdf2a20SBarry Smith } 59ffdf2a20SBarry Smith 607087cfbeSBarry Smith PetscErrorCode SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 619308c137SBarry Smith { 629308c137SBarry Smith PetscErrorCode ierr; 634d4332d5SBarry Smith PetscViewer viewer = (PetscViewer) dummy; 649308c137SBarry Smith const PetscScalar *x,*xl,*xu,*f; 656fd67740SBarry Smith PetscInt i,n,act[2] = {0,0},fact[2],N; 666a9e2d46SJungho Lee /* Number of components that actually hit the bounds (c.f. active variables) */ 676a9e2d46SJungho Lee PetscInt act_bound[2] = {0,0},fact_bound[2]; 68349187a7SBarry Smith PetscReal rnorm,fnorm,zerotolerance = snes->vizerotolerance; 699d1809e2SSatish Balay double tmp; 709308c137SBarry Smith 719308c137SBarry Smith PetscFunctionBegin; 724d4332d5SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 739308c137SBarry Smith ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr); 746fd67740SBarry Smith ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr); 75c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 76c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 779308c137SBarry Smith ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 783f731af5SBarry Smith ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr); 799308c137SBarry Smith 809308c137SBarry Smith rnorm = 0.0; 819308c137SBarry Smith for (i=0; i<n; i++) { 82349187a7SBarry 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]); 83349187a7SBarry Smith else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++; 84349187a7SBarry Smith else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++; 85ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Can never get here"); 869308c137SBarry Smith } 876a9e2d46SJungho Lee 886a9e2d46SJungho Lee for (i=0; i<n; i++) { 89349187a7SBarry Smith if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++; 90349187a7SBarry Smith else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++; 916a9e2d46SJungho Lee } 923f731af5SBarry Smith ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr); 93c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 94c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 959308c137SBarry Smith ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 96b2566f29SBarry Smith ierr = MPIU_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 97b2566f29SBarry Smith ierr = MPIU_Allreduce(act,fact,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 98b2566f29SBarry Smith ierr = MPIU_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 99f137e44dSBarry Smith fnorm = PetscSqrtReal(fnorm); 1006fd67740SBarry Smith 101649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 1029d1809e2SSatish Balay if (snes->ntruebounds) tmp = ((double)(fact[0]+fact[1]))/((double)snes->ntruebounds); 1039d1809e2SSatish Balay else tmp = 0.0; 10462151390SBarry 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); 1056a9e2d46SJungho Lee 106649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 1079308c137SBarry Smith PetscFunctionReturn(0); 1089308c137SBarry Smith } 1099308c137SBarry Smith 1102b4ed282SShri Abhyankar /* 1112b4ed282SShri Abhyankar Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 1122b4ed282SShri 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 1132b4ed282SShri Abhyankar 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 1142b4ed282SShri Abhyankar for this trick. One assumes that the probability that W is in the null space of J is very, very small. 1152b4ed282SShri Abhyankar */ 116ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 1172b4ed282SShri Abhyankar { 1182b4ed282SShri Abhyankar PetscReal a1; 1192b4ed282SShri Abhyankar PetscErrorCode ierr; 120ace3abfcSBarry Smith PetscBool hastranspose; 1212b4ed282SShri Abhyankar 1222b4ed282SShri Abhyankar PetscFunctionBegin; 1232b4ed282SShri Abhyankar *ismin = PETSC_FALSE; 1242b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 1252b4ed282SShri Abhyankar if (hastranspose) { 1262b4ed282SShri Abhyankar /* Compute || J^T F|| */ 1272b4ed282SShri Abhyankar ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 1282b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 1294839bfe8SBarry Smith ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr); 1302b4ed282SShri Abhyankar if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 1312b4ed282SShri Abhyankar } else { 1322b4ed282SShri Abhyankar Vec work; 1332b4ed282SShri Abhyankar PetscScalar result; 1342b4ed282SShri Abhyankar PetscReal wnorm; 1352b4ed282SShri Abhyankar 1360298fd71SBarry Smith ierr = VecSetRandom(W,NULL);CHKERRQ(ierr); 1372b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 1382b4ed282SShri Abhyankar ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 1392b4ed282SShri Abhyankar ierr = MatMult(A,W,work);CHKERRQ(ierr); 1402b4ed282SShri Abhyankar ierr = VecDot(F,work,&result);CHKERRQ(ierr); 1416bf464f9SBarry Smith ierr = VecDestroy(&work);CHKERRQ(ierr); 1422b4ed282SShri Abhyankar a1 = PetscAbsScalar(result)/(fnorm*wnorm); 1434839bfe8SBarry Smith ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr); 1442b4ed282SShri Abhyankar if (a1 < 1.e-4) *ismin = PETSC_TRUE; 1452b4ed282SShri Abhyankar } 1462b4ed282SShri Abhyankar PetscFunctionReturn(0); 1472b4ed282SShri Abhyankar } 1482b4ed282SShri Abhyankar 1492b4ed282SShri Abhyankar /* 1502b4ed282SShri Abhyankar Checks if J^T(F - J*X) = 0 1512b4ed282SShri Abhyankar */ 1522b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 1532b4ed282SShri Abhyankar { 1542b4ed282SShri Abhyankar PetscReal a1,a2; 1552b4ed282SShri Abhyankar PetscErrorCode ierr; 156ace3abfcSBarry Smith PetscBool hastranspose; 1572b4ed282SShri Abhyankar 1582b4ed282SShri Abhyankar PetscFunctionBegin; 1592b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 1602b4ed282SShri Abhyankar if (hastranspose) { 1612b4ed282SShri Abhyankar ierr = MatMult(A,X,W1);CHKERRQ(ierr); 1622b4ed282SShri Abhyankar ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 1632b4ed282SShri Abhyankar 1642b4ed282SShri Abhyankar /* Compute || J^T W|| */ 1652b4ed282SShri Abhyankar ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 1662b4ed282SShri Abhyankar ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 1672b4ed282SShri Abhyankar ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 1682b4ed282SShri Abhyankar if (a1 != 0.0) { 1694839bfe8SBarry Smith ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr); 1702b4ed282SShri Abhyankar } 1712b4ed282SShri Abhyankar } 1722b4ed282SShri Abhyankar PetscFunctionReturn(0); 1732b4ed282SShri Abhyankar } 1742b4ed282SShri Abhyankar 1752b4ed282SShri Abhyankar /* 1768d359177SBarry Smith SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm. 1772b4ed282SShri Abhyankar 1782b4ed282SShri Abhyankar Notes: 1792b4ed282SShri Abhyankar The convergence criterion currently implemented is 1802b4ed282SShri Abhyankar merit < abstol 1812b4ed282SShri Abhyankar merit < rtol*merit_initial 1822b4ed282SShri Abhyankar */ 1838d359177SBarry Smith PetscErrorCode SNESConvergedDefault_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 1842b4ed282SShri Abhyankar { 1852b4ed282SShri Abhyankar PetscErrorCode ierr; 1862b4ed282SShri Abhyankar 1872b4ed282SShri Abhyankar PetscFunctionBegin; 1882b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1892b4ed282SShri Abhyankar PetscValidPointer(reason,6); 1902b4ed282SShri Abhyankar 1912b4ed282SShri Abhyankar *reason = SNES_CONVERGED_ITERATING; 1922b4ed282SShri Abhyankar 1932b4ed282SShri Abhyankar if (!it) { 1942b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 1957fe79bc4SShri Abhyankar snes->ttol = fnorm*snes->rtol; 1962b4ed282SShri Abhyankar } 1977fe79bc4SShri Abhyankar if (fnorm != fnorm) { 1982b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 1992b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FNORM_NAN; 2000d6f27a8SBarry Smith } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) { 2014839bfe8SBarry Smith ierr = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr); 2022b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_ABS; 203e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 2042b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 2052b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FUNCTION_COUNT; 2062b4ed282SShri Abhyankar } 2072b4ed282SShri Abhyankar 2082b4ed282SShri Abhyankar if (it && !*reason) { 2097fe79bc4SShri Abhyankar if (fnorm < snes->ttol) { 2104839bfe8SBarry Smith ierr = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr); 2112b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_RELATIVE; 2122b4ed282SShri Abhyankar } 2132b4ed282SShri Abhyankar } 2142b4ed282SShri Abhyankar PetscFunctionReturn(0); 2152b4ed282SShri Abhyankar } 2162b4ed282SShri Abhyankar 217ee657d29SShri Abhyankar 218c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */ 219c1a5e521SShri Abhyankar /* 220c1a5e521SShri Abhyankar SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n. 221c1a5e521SShri Abhyankar 222c1a5e521SShri Abhyankar Input Parameters: 223c1a5e521SShri Abhyankar . SNES - nonlinear solver context 224c1a5e521SShri Abhyankar 225c1a5e521SShri Abhyankar Output Parameters: 226c1a5e521SShri Abhyankar . X - Bound projected X 227c1a5e521SShri Abhyankar 228c1a5e521SShri Abhyankar */ 229c1a5e521SShri Abhyankar 230c1a5e521SShri Abhyankar PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X) 231c1a5e521SShri Abhyankar { 232c1a5e521SShri Abhyankar PetscErrorCode ierr; 233c1a5e521SShri Abhyankar const PetscScalar *xl,*xu; 234c1a5e521SShri Abhyankar PetscScalar *x; 235c1a5e521SShri Abhyankar PetscInt i,n; 236c1a5e521SShri Abhyankar 237c1a5e521SShri Abhyankar PetscFunctionBegin; 238c1a5e521SShri Abhyankar ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 239c1a5e521SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 240c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 241c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 242c1a5e521SShri Abhyankar 243c1a5e521SShri Abhyankar for (i = 0; i<n; i++) { 24410f5ae6bSBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i]; 24510f5ae6bSBarry Smith else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i]; 246c1a5e521SShri Abhyankar } 247c1a5e521SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 248c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 249c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 250c1a5e521SShri Abhyankar PetscFunctionReturn(0); 251c1a5e521SShri Abhyankar } 252c1a5e521SShri Abhyankar 25369c03620SShri Abhyankar 254693ddf92SShri Abhyankar /* 255693ddf92SShri Abhyankar SNESVIGetActiveSetIndices - Gets the global indices for the active set variables 256693ddf92SShri Abhyankar 257693ddf92SShri Abhyankar Input parameter 258693ddf92SShri Abhyankar . snes - the SNES context 259693ddf92SShri Abhyankar . X - the snes solution vector 260693ddf92SShri Abhyankar . F - the nonlinear function vector 261693ddf92SShri Abhyankar 262693ddf92SShri Abhyankar Output parameter 263693ddf92SShri Abhyankar . ISact - active set index set 264693ddf92SShri Abhyankar */ 265693ddf92SShri Abhyankar PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS *ISact) 266d950fb63SShri Abhyankar { 267d950fb63SShri Abhyankar PetscErrorCode ierr; 268c2fc9fa9SBarry Smith Vec Xl=snes->xl,Xu=snes->xu; 269693ddf92SShri Abhyankar const PetscScalar *x,*f,*xl,*xu; 270693ddf92SShri Abhyankar PetscInt *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0; 271349187a7SBarry Smith PetscReal zerotolerance = snes->vizerotolerance; 272d950fb63SShri Abhyankar 273d950fb63SShri Abhyankar PetscFunctionBegin; 274d950fb63SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 275d950fb63SShri Abhyankar ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr); 276fe835674SShri Abhyankar ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 277fe835674SShri Abhyankar ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr); 278fe835674SShri Abhyankar ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr); 279fe835674SShri Abhyankar ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 280693ddf92SShri Abhyankar /* Compute active set size */ 281d950fb63SShri Abhyankar for (i=0; i < nlocal;i++) { 282349187a7SBarry 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++; 283d950fb63SShri Abhyankar } 284693ddf92SShri Abhyankar 285785e854fSJed Brown ierr = PetscMalloc1(nloc_isact,&idx_act);CHKERRQ(ierr); 286d950fb63SShri Abhyankar 287693ddf92SShri Abhyankar /* Set active set indices */ 288d950fb63SShri Abhyankar for (i=0; i < nlocal; i++) { 289349187a7SBarry 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; 290d950fb63SShri Abhyankar } 291d950fb63SShri Abhyankar 292693ddf92SShri Abhyankar /* Create active set IS */ 293ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes),nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr); 294d950fb63SShri Abhyankar 295fe835674SShri Abhyankar ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 296fe835674SShri Abhyankar ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr); 297fe835674SShri Abhyankar ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr); 298fe835674SShri Abhyankar ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 299d950fb63SShri Abhyankar PetscFunctionReturn(0); 300d950fb63SShri Abhyankar } 301d950fb63SShri Abhyankar 302693ddf92SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS *ISact,IS *ISinact) 303693ddf92SShri Abhyankar { 304693ddf92SShri Abhyankar PetscErrorCode ierr; 305077aedafSJed Brown PetscInt rstart,rend; 306693ddf92SShri Abhyankar 307693ddf92SShri Abhyankar PetscFunctionBegin; 308693ddf92SShri Abhyankar ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr); 309077aedafSJed Brown ierr = VecGetOwnershipRange(X,&rstart,&rend);CHKERRQ(ierr); 310077aedafSJed Brown ierr = ISComplement(*ISact,rstart,rend,ISinact);CHKERRQ(ierr); 311693ddf92SShri Abhyankar PetscFunctionReturn(0); 312693ddf92SShri Abhyankar } 313693ddf92SShri Abhyankar 31410f5ae6bSBarry Smith PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm) 315fe835674SShri Abhyankar { 316fe835674SShri Abhyankar PetscErrorCode ierr; 317fe835674SShri Abhyankar const PetscScalar *x,*xl,*xu,*f; 318fe835674SShri Abhyankar PetscInt i,n; 319349187a7SBarry Smith PetscReal rnorm,zerotolerance = snes->vizerotolerance;; 320fe835674SShri Abhyankar 321fe835674SShri Abhyankar PetscFunctionBegin; 322fe835674SShri Abhyankar ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 323c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 324c2fc9fa9SBarry Smith ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 325fe835674SShri Abhyankar ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 326fe835674SShri Abhyankar ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 327fe835674SShri Abhyankar rnorm = 0.0; 328fe835674SShri Abhyankar for (i=0; i<n; i++) { 329349187a7SBarry 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]); 3308f918527SKarl Rupp } 331fe835674SShri Abhyankar ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 332c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 333c2fc9fa9SBarry Smith ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 334fe835674SShri Abhyankar ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 335b2566f29SBarry Smith ierr = MPIU_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 33662d1f40fSBarry Smith *fnorm = PetscSqrtReal(*fnorm); 337fe835674SShri Abhyankar PetscFunctionReturn(0); 338fe835674SShri Abhyankar } 339fe835674SShri Abhyankar 34008da532bSDmitry Karpeev PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu) 34108da532bSDmitry Karpeev { 34208da532bSDmitry Karpeev PetscErrorCode ierr; 3436e111a19SKarl Rupp 34408da532bSDmitry Karpeev PetscFunctionBegin; 34508da532bSDmitry Karpeev ierr = DMComputeVariableBounds(snes->dm, xl, xu);CHKERRQ(ierr); 34608da532bSDmitry Karpeev PetscFunctionReturn(0); 34708da532bSDmitry Karpeev } 34808da532bSDmitry Karpeev 3492f63a38cSShri Abhyankar 3502b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 3512b4ed282SShri Abhyankar /* 352c2fc9fa9SBarry Smith SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up 3532b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 3542b4ed282SShri Abhyankar 3552b4ed282SShri Abhyankar Input Parameter: 3562b4ed282SShri Abhyankar . snes - the SNES context 3572b4ed282SShri Abhyankar 3582b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 3592b4ed282SShri Abhyankar 3602b4ed282SShri Abhyankar Notes: 3612b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 3622b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 3632b4ed282SShri Abhyankar the call to SNESSolve(). 3642b4ed282SShri Abhyankar */ 3652b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes) 3662b4ed282SShri Abhyankar { 3672b4ed282SShri Abhyankar PetscErrorCode ierr; 3682b4ed282SShri Abhyankar PetscInt i_start[3],i_end[3]; 3692b4ed282SShri Abhyankar 3702b4ed282SShri Abhyankar PetscFunctionBegin; 371fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,1);CHKERRQ(ierr); 3726cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 3732b4ed282SShri Abhyankar 37408da532bSDmitry Karpeev if (!snes->ops->computevariablebounds && snes->dm) { 375a201590fSDmitry Karpeev PetscBool flag; 376a201590fSDmitry Karpeev ierr = DMHasVariableBounds(snes->dm, &flag);CHKERRQ(ierr); 37736b2a9f3SBarry Smith if (flag) { 37808da532bSDmitry Karpeev snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds; 37908da532bSDmitry Karpeev } 38036b2a9f3SBarry Smith } 381a201590fSDmitry Karpeev if (!snes->usersetbounds) { 382c2fc9fa9SBarry Smith if (snes->ops->computevariablebounds) { 383c2fc9fa9SBarry Smith if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);} 384c2fc9fa9SBarry Smith if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);} 385c2fc9fa9SBarry Smith ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr); 3861aa26658SKarl Rupp } else if (!snes->xl && !snes->xu) { 3872176524cSBarry Smith /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */ 388c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr); 389e270355aSBarry Smith ierr = VecSet(snes->xl,PETSC_NINFINITY);CHKERRQ(ierr); 390c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr); 391e270355aSBarry Smith ierr = VecSet(snes->xu,PETSC_INFINITY);CHKERRQ(ierr); 3922b4ed282SShri Abhyankar } else { 3932b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 3942b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 395c2fc9fa9SBarry Smith ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr); 396c2fc9fa9SBarry Smith ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr); 3972b4ed282SShri 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])) 3982b4ed282SShri 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."); 3992b4ed282SShri Abhyankar } 400a201590fSDmitry Karpeev } 4012b4ed282SShri Abhyankar PetscFunctionReturn(0); 4022b4ed282SShri Abhyankar } 4032b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 4042176524cSBarry Smith PetscErrorCode SNESReset_VI(SNES snes) 4052176524cSBarry Smith { 4062176524cSBarry Smith PetscErrorCode ierr; 4072176524cSBarry Smith 4082176524cSBarry Smith PetscFunctionBegin; 409c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xl);CHKERRQ(ierr); 410c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xu);CHKERRQ(ierr); 4112d6615e8SDmitry Karpeev snes->usersetbounds = PETSC_FALSE; 4122176524cSBarry Smith PetscFunctionReturn(0); 4132176524cSBarry Smith } 4142176524cSBarry Smith 4152b4ed282SShri Abhyankar /* 4162b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 4172b4ed282SShri Abhyankar with SNESCreate_VI(). 4182b4ed282SShri Abhyankar 4192b4ed282SShri Abhyankar Input Parameter: 4202b4ed282SShri Abhyankar . snes - the SNES context 4212b4ed282SShri Abhyankar 4222b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 4232b4ed282SShri Abhyankar */ 4242b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes) 4252b4ed282SShri Abhyankar { 4262b4ed282SShri Abhyankar PetscErrorCode ierr; 4272b4ed282SShri Abhyankar 4282b4ed282SShri Abhyankar PetscFunctionBegin; 4292b4ed282SShri Abhyankar ierr = PetscFree(snes->data);CHKERRQ(ierr); 4302b4ed282SShri Abhyankar 4312b4ed282SShri Abhyankar /* clear composed functions */ 432bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSet_C",NULL);CHKERRQ(ierr); 433dcf2fd19SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESLineSearchSetDefaultMonitor_C",NULL);CHKERRQ(ierr); 4342b4ed282SShri Abhyankar PetscFunctionReturn(0); 4352b4ed282SShri Abhyankar } 4367fe79bc4SShri Abhyankar 4375559b345SBarry Smith /*@ 4382b4ed282SShri Abhyankar SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 4392b4ed282SShri Abhyankar 4402b4ed282SShri Abhyankar Input Parameters: 441*a2b725a8SWilliam Gropp + snes - the SNES context. 4422b4ed282SShri Abhyankar . xl - lower bound. 443*a2b725a8SWilliam Gropp - xu - upper bound. 4442b4ed282SShri Abhyankar 4452b4ed282SShri Abhyankar Notes: 4462b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 44729eed3a4SBarry Smith PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp(). 44884c105d7SBarry Smith 4492bd2b0e6SSatish Balay Level: advanced 4502bd2b0e6SSatish Balay 4515559b345SBarry Smith @*/ 45269c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 4532b4ed282SShri Abhyankar { 45461589011SJed Brown PetscErrorCode ierr,(*f)(SNES,Vec,Vec); 4552b4ed282SShri Abhyankar 4562b4ed282SShri Abhyankar PetscFunctionBegin; 4572b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4582b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 4592b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 4600005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",&f);CHKERRQ(ierr); 46156e5e3feSPatrick Farrell if (!f) { 46256e5e3feSPatrick Farrell ierr = SNESVISetVariableBounds_VI(snes, xl, xu);CHKERRQ(ierr); 46356e5e3feSPatrick Farrell } else { 46461589011SJed Brown ierr = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr); 46556e5e3feSPatrick Farrell } 466a201590fSDmitry Karpeev snes->usersetbounds = PETSC_TRUE; 46761589011SJed Brown PetscFunctionReturn(0); 46861589011SJed Brown } 46961589011SJed Brown 47061589011SJed Brown PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu) 47161589011SJed Brown { 47261589011SJed Brown PetscErrorCode ierr; 47361589011SJed Brown const PetscScalar *xxl,*xxu; 47461589011SJed Brown PetscInt i,n, cnt = 0; 47561589011SJed Brown 47661589011SJed Brown PetscFunctionBegin; 4770298fd71SBarry Smith ierr = SNESGetFunction(snes,&snes->vec_func,NULL,NULL);CHKERRQ(ierr); 478a63bb30eSJed Brown if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first"); 479077aedafSJed Brown { 480077aedafSJed Brown PetscInt xlN,xuN,N; 481077aedafSJed Brown ierr = VecGetSize(xl,&xlN);CHKERRQ(ierr); 482077aedafSJed Brown ierr = VecGetSize(xu,&xuN);CHKERRQ(ierr); 483077aedafSJed Brown ierr = VecGetSize(snes->vec_func,&N);CHKERRQ(ierr); 484077aedafSJed Brown if (xlN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xlN,N); 485077aedafSJed Brown if (xuN != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xuN,N); 486077aedafSJed Brown } 4872176524cSBarry Smith ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr); 4882176524cSBarry Smith ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr); 489c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xl);CHKERRQ(ierr); 490c2fc9fa9SBarry Smith ierr = VecDestroy(&snes->xu);CHKERRQ(ierr); 491c2fc9fa9SBarry Smith snes->xl = xl; 492c2fc9fa9SBarry Smith snes->xu = xu; 4936fd67740SBarry Smith ierr = VecGetLocalSize(xl,&n);CHKERRQ(ierr); 4946fd67740SBarry Smith ierr = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr); 4956fd67740SBarry Smith ierr = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr); 496e270355aSBarry Smith for (i=0; i<n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY)); 4971aa26658SKarl Rupp 498b2566f29SBarry Smith ierr = MPIU_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 4996fd67740SBarry Smith ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr); 5006fd67740SBarry Smith ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr); 5012b4ed282SShri Abhyankar PetscFunctionReturn(0); 5022b4ed282SShri Abhyankar } 50392c02d66SPeter Brune 5044416b707SBarry Smith PetscErrorCode SNESSetFromOptions_VI(PetscOptionItems *PetscOptionsObject,SNES snes) 5052b4ed282SShri Abhyankar { 5062b4ed282SShri Abhyankar PetscErrorCode ierr; 5078afaa268SBarry Smith PetscBool flg = PETSC_FALSE; 508639ea3faSPeter Brune SNESLineSearch linesearch; 5092b4ed282SShri Abhyankar 5102b4ed282SShri Abhyankar PetscFunctionBegin; 511e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES VI options");CHKERRQ(ierr); 512349187a7SBarry 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); 513ffdf2a20SBarry Smith ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","SNESMonitorResidual",flg,&flg,NULL);CHKERRQ(ierr); 514c2fc9fa9SBarry Smith if (flg) { 5154d4332d5SBarry Smith ierr = SNESMonitorSet(snes,SNESMonitorVI,PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),NULL);CHKERRQ(ierr); 5162b4ed282SShri Abhyankar } 517de34d3e9SBarry Smith flg = PETSC_FALSE; 518ffdf2a20SBarry Smith ierr = PetscOptionsBool("-snes_vi_monitor_residual","Monitor residual all non-active variables; using zero for active constraints","SNESMonitorVIResidual",flg,&flg,NULL);CHKERRQ(ierr); 519ffdf2a20SBarry Smith if (flg) { 520ffdf2a20SBarry Smith ierr = SNESMonitorSet(snes,SNESVIMonitorResidual,PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)),NULL);CHKERRQ(ierr); 521ffdf2a20SBarry Smith } 522639ea3faSPeter Brune if (!snes->linesearch) { 5237601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 524639ea3faSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr); 525639ea3faSPeter Brune ierr = SNESLineSearchBTSetAlpha(linesearch, 0.0);CHKERRQ(ierr); 526639ea3faSPeter Brune } 527c2fc9fa9SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 528ed70e96aSJungho Lee PetscFunctionReturn(0); 529ed70e96aSJungho Lee } 530