1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 21e25c274SJed Brown #include <petscdm.h> 3d0af7cd3SBarry Smith 42176524cSBarry Smith /*@C 5f0b84518SBarry Smith SNESVISetComputeVariableBounds - Sets a function that is called to compute the variable bounds. This allows solving 6f0b84518SBarry Smith (differential) variable inequalities. For example, restricting pressure to be non-negative. 72176524cSBarry Smith 87a7aea1fSJed Brown Input parameter: 92176524cSBarry Smith + snes - the SNES context 102176524cSBarry Smith - compute - computes the bounds 112176524cSBarry Smith 122bd2b0e6SSatish Balay Level: advanced 132176524cSBarry Smith 14f0b84518SBarry Smith Note: 15f0b84518SBarry Smith Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers. 16f0b84518SBarry Smith 17f0b84518SBarry Smith .seealso: `SNESVISetVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, 18f0b84518SBarry Smith 'SNESSetType()` 19c1c3a0ecSBarry Smith 20aab9d709SJed Brown @*/ 219371c9d4SSatish Balay PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES, Vec, Vec)) { 225f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(SNES, PetscErrorCode(*)(SNES, Vec, Vec)); 232176524cSBarry Smith 242176524cSBarry Smith PetscFunctionBegin; 2561589011SJed Brown PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 269566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", &f)); 27cac4c232SBarry Smith if (f) PetscUseMethod(snes, "SNESVISetComputeVariableBounds_C", (SNES, PetscErrorCode(*)(SNES, Vec, Vec)), (snes, compute)); 289566063dSJacob Faibussowitsch else PetscCall(SNESVISetComputeVariableBounds_VI(snes, compute)); 292176524cSBarry Smith PetscFunctionReturn(0); 302176524cSBarry Smith } 312176524cSBarry Smith 329371c9d4SSatish Balay PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes, SNESVIComputeVariableBoundsFunction compute) { 3361589011SJed Brown PetscFunctionBegin; 3461589011SJed Brown snes->ops->computevariablebounds = compute; 3561589011SJed Brown PetscFunctionReturn(0); 3661589011SJed Brown } 372176524cSBarry Smith 383c0c59f3SBarry Smith /* --------------------------------------------------------------------------------------------------------*/ 392b4ed282SShri Abhyankar 409371c9d4SSatish Balay PetscErrorCode SNESVIMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy) { 41ffdf2a20SBarry Smith Vec X, F, Finactive; 42ffdf2a20SBarry Smith IS isactive; 43ffdf2a20SBarry Smith PetscViewer viewer = (PetscViewer)dummy; 44ffdf2a20SBarry Smith 45ffdf2a20SBarry Smith PetscFunctionBegin; 469566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 479566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &X)); 489566063dSJacob Faibussowitsch PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive)); 499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &Finactive)); 509566063dSJacob Faibussowitsch PetscCall(VecCopy(F, Finactive)); 519566063dSJacob Faibussowitsch PetscCall(VecISSet(Finactive, isactive, 0.0)); 529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isactive)); 539566063dSJacob Faibussowitsch PetscCall(VecView(Finactive, viewer)); 549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Finactive)); 55ffdf2a20SBarry Smith PetscFunctionReturn(0); 56ffdf2a20SBarry Smith } 57ffdf2a20SBarry Smith 589371c9d4SSatish Balay PetscErrorCode SNESMonitorVI(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy) { 594d4332d5SBarry Smith PetscViewer viewer = (PetscViewer)dummy; 609308c137SBarry Smith const PetscScalar *x, *xl, *xu, *f; 616fd67740SBarry Smith PetscInt i, n, act[2] = {0, 0}, fact[2], N; 626a9e2d46SJungho Lee /* Number of components that actually hit the bounds (c.f. active variables) */ 636a9e2d46SJungho Lee PetscInt act_bound[2] = {0, 0}, fact_bound[2]; 64349187a7SBarry Smith PetscReal rnorm, fnorm, zerotolerance = snes->vizerotolerance; 659d1809e2SSatish Balay double tmp; 669308c137SBarry Smith 679308c137SBarry Smith PetscFunctionBegin; 684d4332d5SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4); 699566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 709566063dSJacob Faibussowitsch PetscCall(VecGetSize(snes->vec_sol, &N)); 719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xl, &xl)); 729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xu, &xu)); 739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->vec_sol, &x)); 749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->vec_func, &f)); 759308c137SBarry Smith 769308c137SBarry Smith rnorm = 0.0; 779308c137SBarry Smith for (i = 0; i < n; i++) { 78349187a7SBarry 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]); 79349187a7SBarry Smith else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++; 80349187a7SBarry Smith else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++; 81ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Can never get here"); 829308c137SBarry Smith } 836a9e2d46SJungho Lee 846a9e2d46SJungho Lee for (i = 0; i < n; i++) { 85349187a7SBarry Smith if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++; 86349187a7SBarry Smith else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++; 876a9e2d46SJungho Lee } 889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->vec_func, &f)); 899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->vec_sol, &x)); 921c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&rnorm, &fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes))); 931c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(act, fact, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 941c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(act_bound, fact_bound, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 95f137e44dSBarry Smith fnorm = PetscSqrtReal(fnorm); 966fd67740SBarry Smith 979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel)); 989d1809e2SSatish Balay if (snes->ntruebounds) tmp = ((double)(fact[0] + fact[1])) / ((double)snes->ntruebounds); 999d1809e2SSatish Balay else tmp = 0.0; 10063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES VI Function norm %g Active lower constraints %" PetscInt_FMT "/%" PetscInt_FMT " upper constraints %" PetscInt_FMT "/%" PetscInt_FMT " 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)); 1016a9e2d46SJungho Lee 1029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel)); 1039308c137SBarry Smith PetscFunctionReturn(0); 1049308c137SBarry Smith } 1059308c137SBarry Smith 1062b4ed282SShri Abhyankar /* 1072b4ed282SShri Abhyankar Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 1082b4ed282SShri 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 1092b4ed282SShri Abhyankar 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 1102b4ed282SShri Abhyankar for this trick. One assumes that the probability that W is in the null space of J is very, very small. 1112b4ed282SShri Abhyankar */ 1129371c9d4SSatish Balay PetscErrorCode SNESVICheckLocalMin_Private(SNES snes, Mat A, Vec F, Vec W, PetscReal fnorm, PetscBool *ismin) { 1132b4ed282SShri Abhyankar PetscReal a1; 114ace3abfcSBarry Smith PetscBool hastranspose; 1152b4ed282SShri Abhyankar 1162b4ed282SShri Abhyankar PetscFunctionBegin; 1172b4ed282SShri Abhyankar *ismin = PETSC_FALSE; 1189566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 1192b4ed282SShri Abhyankar if (hastranspose) { 1202b4ed282SShri Abhyankar /* Compute || J^T F|| */ 1219566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, F, W)); 1229566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &a1)); 1239566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "|| J^T F|| %g near zero implies found a local minimum\n", (double)(a1 / fnorm))); 1242b4ed282SShri Abhyankar if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE; 1252b4ed282SShri Abhyankar } else { 1262b4ed282SShri Abhyankar Vec work; 1272b4ed282SShri Abhyankar PetscScalar result; 1282b4ed282SShri Abhyankar PetscReal wnorm; 1292b4ed282SShri Abhyankar 1309566063dSJacob Faibussowitsch PetscCall(VecSetRandom(W, NULL)); 1319566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 1329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(W, &work)); 1339566063dSJacob Faibussowitsch PetscCall(MatMult(A, W, work)); 1349566063dSJacob Faibussowitsch PetscCall(VecDot(F, work, &result)); 1359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 1362b4ed282SShri Abhyankar a1 = PetscAbsScalar(result) / (fnorm * wnorm); 1379566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n", (double)a1)); 1382b4ed282SShri Abhyankar if (a1 < 1.e-4) *ismin = PETSC_TRUE; 1392b4ed282SShri Abhyankar } 1402b4ed282SShri Abhyankar PetscFunctionReturn(0); 1412b4ed282SShri Abhyankar } 1422b4ed282SShri Abhyankar 1432b4ed282SShri Abhyankar /* 1442b4ed282SShri Abhyankar Checks if J^T(F - J*X) = 0 1452b4ed282SShri Abhyankar */ 1469371c9d4SSatish Balay PetscErrorCode SNESVICheckResidual_Private(SNES snes, Mat A, Vec F, Vec X, Vec W1, Vec W2) { 1472b4ed282SShri Abhyankar PetscReal a1, a2; 148ace3abfcSBarry Smith PetscBool hastranspose; 1492b4ed282SShri Abhyankar 1502b4ed282SShri Abhyankar PetscFunctionBegin; 1519566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 1522b4ed282SShri Abhyankar if (hastranspose) { 1539566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, W1)); 1549566063dSJacob Faibussowitsch PetscCall(VecAXPY(W1, -1.0, F)); 1552b4ed282SShri Abhyankar 1562b4ed282SShri Abhyankar /* Compute || J^T W|| */ 1579566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, W1, W2)); 1589566063dSJacob Faibussowitsch PetscCall(VecNorm(W1, NORM_2, &a1)); 1599566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &a2)); 160*48a46eb9SPierre Jolivet if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n", (double)(a2 / a1))); 1612b4ed282SShri Abhyankar } 1622b4ed282SShri Abhyankar PetscFunctionReturn(0); 1632b4ed282SShri Abhyankar } 1642b4ed282SShri Abhyankar 1652b4ed282SShri Abhyankar /* 1668d359177SBarry Smith SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm. 1672b4ed282SShri Abhyankar 1682b4ed282SShri Abhyankar Notes: 1692b4ed282SShri Abhyankar The convergence criterion currently implemented is 1702b4ed282SShri Abhyankar merit < abstol 1712b4ed282SShri Abhyankar merit < rtol*merit_initial 1722b4ed282SShri Abhyankar */ 1739371c9d4SSatish Balay PetscErrorCode SNESConvergedDefault_VI(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gradnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) { 1742b4ed282SShri Abhyankar PetscFunctionBegin; 1752b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1762b4ed282SShri Abhyankar PetscValidPointer(reason, 6); 1772b4ed282SShri Abhyankar 1782b4ed282SShri Abhyankar *reason = SNES_CONVERGED_ITERATING; 1792b4ed282SShri Abhyankar 1802b4ed282SShri Abhyankar if (!it) { 1812b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 1827fe79bc4SShri Abhyankar snes->ttol = fnorm * snes->rtol; 1832b4ed282SShri Abhyankar } 1847fe79bc4SShri Abhyankar if (fnorm != fnorm) { 1859566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n")); 1862b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FNORM_NAN; 1870d6f27a8SBarry Smith } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) { 1889566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g\n", (double)fnorm, (double)snes->abstol)); 1892b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_ABS; 190e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 19163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs)); 1922b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FUNCTION_COUNT; 1932b4ed282SShri Abhyankar } 1942b4ed282SShri Abhyankar 1952b4ed282SShri Abhyankar if (it && !*reason) { 1967fe79bc4SShri Abhyankar if (fnorm < snes->ttol) { 1979566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g (relative tolerance)\n", (double)fnorm, (double)snes->ttol)); 1982b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_RELATIVE; 1992b4ed282SShri Abhyankar } 2002b4ed282SShri Abhyankar } 2012b4ed282SShri Abhyankar PetscFunctionReturn(0); 2022b4ed282SShri Abhyankar } 2032b4ed282SShri Abhyankar 204c1a5e521SShri Abhyankar /* -------------------------------------------------------------------------- */ 205c1a5e521SShri Abhyankar /* 206c1a5e521SShri Abhyankar SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n. 207c1a5e521SShri Abhyankar 208c1a5e521SShri Abhyankar Input Parameters: 209c1a5e521SShri Abhyankar . SNES - nonlinear solver context 210c1a5e521SShri Abhyankar 211c1a5e521SShri Abhyankar Output Parameters: 212c1a5e521SShri Abhyankar . X - Bound projected X 213c1a5e521SShri Abhyankar 214c1a5e521SShri Abhyankar */ 215c1a5e521SShri Abhyankar 2169371c9d4SSatish Balay PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X) { 217c1a5e521SShri Abhyankar const PetscScalar *xl, *xu; 218c1a5e521SShri Abhyankar PetscScalar *x; 219c1a5e521SShri Abhyankar PetscInt i, n; 220c1a5e521SShri Abhyankar 221c1a5e521SShri Abhyankar PetscFunctionBegin; 2229566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 2239566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 2249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xl, &xl)); 2259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xu, &xu)); 226c1a5e521SShri Abhyankar 227c1a5e521SShri Abhyankar for (i = 0; i < n; i++) { 22810f5ae6bSBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i]; 22910f5ae6bSBarry Smith else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i]; 230c1a5e521SShri Abhyankar } 2319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 2329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 2339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 234c1a5e521SShri Abhyankar PetscFunctionReturn(0); 235c1a5e521SShri Abhyankar } 236c1a5e521SShri Abhyankar 237693ddf92SShri Abhyankar /* 238693ddf92SShri Abhyankar SNESVIGetActiveSetIndices - Gets the global indices for the active set variables 239693ddf92SShri Abhyankar 2407a7aea1fSJed Brown Input parameter: 241693ddf92SShri Abhyankar . snes - the SNES context 242693ddf92SShri Abhyankar . X - the snes solution vector 243693ddf92SShri Abhyankar . F - the nonlinear function vector 244693ddf92SShri Abhyankar 2457a7aea1fSJed Brown Output parameter: 246693ddf92SShri Abhyankar . ISact - active set index set 247693ddf92SShri Abhyankar */ 2489371c9d4SSatish Balay PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact) { 249c2fc9fa9SBarry Smith Vec Xl = snes->xl, Xu = snes->xu; 250693ddf92SShri Abhyankar const PetscScalar *x, *f, *xl, *xu; 251693ddf92SShri Abhyankar PetscInt *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0; 252349187a7SBarry Smith PetscReal zerotolerance = snes->vizerotolerance; 253d950fb63SShri Abhyankar 254d950fb63SShri Abhyankar PetscFunctionBegin; 2559566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nlocal)); 2569566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh)); 2579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xl, &xl)); 2599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xu, &xu)); 2609566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &f)); 261693ddf92SShri Abhyankar /* Compute active set size */ 262d950fb63SShri Abhyankar for (i = 0; i < nlocal; i++) { 263349187a7SBarry 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++; 264d950fb63SShri Abhyankar } 265693ddf92SShri Abhyankar 2669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc_isact, &idx_act)); 267d950fb63SShri Abhyankar 268693ddf92SShri Abhyankar /* Set active set indices */ 269d950fb63SShri Abhyankar for (i = 0; i < nlocal; i++) { 270349187a7SBarry 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; 271d950fb63SShri Abhyankar } 272d950fb63SShri Abhyankar 273693ddf92SShri Abhyankar /* Create active set IS */ 2749566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact)); 275d950fb63SShri Abhyankar 2769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xl, &xl)); 2789566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xu, &xu)); 2799566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &f)); 280d950fb63SShri Abhyankar PetscFunctionReturn(0); 281d950fb63SShri Abhyankar } 282d950fb63SShri Abhyankar 2839371c9d4SSatish Balay PetscErrorCode SNESVICreateIndexSets_RS(SNES snes, Vec X, Vec F, IS *ISact, IS *ISinact) { 284077aedafSJed Brown PetscInt rstart, rend; 285693ddf92SShri Abhyankar 286693ddf92SShri Abhyankar PetscFunctionBegin; 2879566063dSJacob Faibussowitsch PetscCall(SNESVIGetActiveSetIS(snes, X, F, ISact)); 2889566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &rstart, &rend)); 2899566063dSJacob Faibussowitsch PetscCall(ISComplement(*ISact, rstart, rend, ISinact)); 290693ddf92SShri Abhyankar PetscFunctionReturn(0); 291693ddf92SShri Abhyankar } 292693ddf92SShri Abhyankar 2939371c9d4SSatish Balay PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm) { 294fe835674SShri Abhyankar const PetscScalar *x, *xl, *xu, *f; 295fe835674SShri Abhyankar PetscInt i, n; 296feb237baSPierre Jolivet PetscReal rnorm, zerotolerance = snes->vizerotolerance; 297fe835674SShri Abhyankar 298fe835674SShri Abhyankar PetscFunctionBegin; 2999566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 3009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xl, &xl)); 3019566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xu, &xu)); 3029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 3039566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &f)); 304fe835674SShri Abhyankar rnorm = 0.0; 305fe835674SShri Abhyankar for (i = 0; i < n; i++) { 306349187a7SBarry 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]); 3078f918527SKarl Rupp } 3089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &f)); 3099566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 3109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 3119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 3121c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes))); 31362d1f40fSBarry Smith *fnorm = PetscSqrtReal(*fnorm); 314fe835674SShri Abhyankar PetscFunctionReturn(0); 315fe835674SShri Abhyankar } 316fe835674SShri Abhyankar 3179371c9d4SSatish Balay PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu) { 31808da532bSDmitry Karpeev PetscFunctionBegin; 3199566063dSJacob Faibussowitsch PetscCall(DMComputeVariableBounds(snes->dm, xl, xu)); 32008da532bSDmitry Karpeev PetscFunctionReturn(0); 32108da532bSDmitry Karpeev } 32208da532bSDmitry Karpeev 3232b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 3242b4ed282SShri Abhyankar /* 325c2fc9fa9SBarry Smith SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up 3262b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 3272b4ed282SShri Abhyankar 3282b4ed282SShri Abhyankar Input Parameter: 3292b4ed282SShri Abhyankar . snes - the SNES context 3302b4ed282SShri Abhyankar 3312b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 3322b4ed282SShri Abhyankar 3332b4ed282SShri Abhyankar Notes: 3342b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 3352b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 3362b4ed282SShri Abhyankar the call to SNESSolve(). 3372b4ed282SShri Abhyankar */ 3389371c9d4SSatish Balay PetscErrorCode SNESSetUp_VI(SNES snes) { 3392b4ed282SShri Abhyankar PetscInt i_start[3], i_end[3]; 3402b4ed282SShri Abhyankar 3412b4ed282SShri Abhyankar PetscFunctionBegin; 3429566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 1)); 3439566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 3442b4ed282SShri Abhyankar 34508da532bSDmitry Karpeev if (!snes->ops->computevariablebounds && snes->dm) { 346a201590fSDmitry Karpeev PetscBool flag; 3479566063dSJacob Faibussowitsch PetscCall(DMHasVariableBounds(snes->dm, &flag)); 3489371c9d4SSatish Balay if (flag) { snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds; } 34936b2a9f3SBarry Smith } 350a201590fSDmitry Karpeev if (!snes->usersetbounds) { 351c2fc9fa9SBarry Smith if (snes->ops->computevariablebounds) { 3529566063dSJacob Faibussowitsch if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl)); 3539566063dSJacob Faibussowitsch if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu)); 354dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu); 3551aa26658SKarl Rupp } else if (!snes->xl && !snes->xu) { 3562176524cSBarry Smith /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */ 3579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &snes->xl)); 3589566063dSJacob Faibussowitsch PetscCall(VecSet(snes->xl, PETSC_NINFINITY)); 3599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &snes->xu)); 3609566063dSJacob Faibussowitsch PetscCall(VecSet(snes->xu, PETSC_INFINITY)); 3612b4ed282SShri Abhyankar } else { 3622b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 3639566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end)); 3649566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1)); 3659566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2)); 3662b4ed282SShri 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])) 3672b4ed282SShri 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."); 3682b4ed282SShri Abhyankar } 369a201590fSDmitry Karpeev } 3702b4ed282SShri Abhyankar PetscFunctionReturn(0); 3712b4ed282SShri Abhyankar } 3722b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 3739371c9d4SSatish Balay PetscErrorCode SNESReset_VI(SNES snes) { 3742176524cSBarry Smith PetscFunctionBegin; 3759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xl)); 3769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xu)); 3772d6615e8SDmitry Karpeev snes->usersetbounds = PETSC_FALSE; 3782176524cSBarry Smith PetscFunctionReturn(0); 3792176524cSBarry Smith } 3802176524cSBarry Smith 3812b4ed282SShri Abhyankar /* 3822b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 3832b4ed282SShri Abhyankar with SNESCreate_VI(). 3842b4ed282SShri Abhyankar 3852b4ed282SShri Abhyankar Input Parameter: 3862b4ed282SShri Abhyankar . snes - the SNES context 3872b4ed282SShri Abhyankar 3882b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 3892b4ed282SShri Abhyankar */ 3909371c9d4SSatish Balay PetscErrorCode SNESDestroy_VI(SNES snes) { 3912b4ed282SShri Abhyankar PetscFunctionBegin; 3929566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3932b4ed282SShri Abhyankar 3942b4ed282SShri Abhyankar /* clear composed functions */ 3952e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL)); 3962e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL)); 3972b4ed282SShri Abhyankar PetscFunctionReturn(0); 3982b4ed282SShri Abhyankar } 3997fe79bc4SShri Abhyankar 4005559b345SBarry Smith /*@ 401f0b84518SBarry Smith SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving 402f0b84518SBarry Smith (differential) variable inequalities. For example, restricting pressure to be non-negative. 4032b4ed282SShri Abhyankar 4042b4ed282SShri Abhyankar Input Parameters: 405a2b725a8SWilliam Gropp + snes - the SNES context. 4062b4ed282SShri Abhyankar . xl - lower bound. 407a2b725a8SWilliam Gropp - xu - upper bound. 4082b4ed282SShri Abhyankar 4092b4ed282SShri Abhyankar Notes: 4102b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 41129eed3a4SBarry Smith PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp(). 41284c105d7SBarry Smith 413f0b84518SBarry Smith Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers. 414f0b84518SBarry Smith 4152bd2b0e6SSatish Balay Level: advanced 4162bd2b0e6SSatish Balay 417f0b84518SBarry Smith .seealso: `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, SNESVINEWTONRSLS, SNESVINEWTONSSLS, 'SNESSetType()` 4185559b345SBarry Smith @*/ 4199371c9d4SSatish Balay PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) { 4205f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(SNES, Vec, Vec); 4212b4ed282SShri Abhyankar 4222b4ed282SShri Abhyankar PetscFunctionBegin; 4232b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 4242b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl, VEC_CLASSID, 2); 4252b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu, VEC_CLASSID, 3); 4269566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f)); 427cac4c232SBarry Smith if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu)); 4289566063dSJacob Faibussowitsch else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu)); 429a201590fSDmitry Karpeev snes->usersetbounds = PETSC_TRUE; 43061589011SJed Brown PetscFunctionReturn(0); 43161589011SJed Brown } 43261589011SJed Brown 4339371c9d4SSatish Balay PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu) { 43461589011SJed Brown const PetscScalar *xxl, *xxu; 43561589011SJed Brown PetscInt i, n, cnt = 0; 43661589011SJed Brown 43761589011SJed Brown PetscFunctionBegin; 4389566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL)); 4395f80ce2aSJacob Faibussowitsch PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first"); 440077aedafSJed Brown { 441077aedafSJed Brown PetscInt xlN, xuN, N; 4429566063dSJacob Faibussowitsch PetscCall(VecGetSize(xl, &xlN)); 4439566063dSJacob Faibussowitsch PetscCall(VecGetSize(xu, &xuN)); 4449566063dSJacob Faibussowitsch PetscCall(VecGetSize(snes->vec_func, &N)); 44563a3b9bcSJacob Faibussowitsch PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N); 44663a3b9bcSJacob Faibussowitsch PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N); 447077aedafSJed Brown } 4489566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xl)); 4499566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xu)); 4509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xl)); 4519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xu)); 452c2fc9fa9SBarry Smith snes->xl = xl; 453c2fc9fa9SBarry Smith snes->xu = xu; 4549566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xl, &n)); 4559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xl, &xxl)); 4569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xu, &xxu)); 457e270355aSBarry Smith for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY)); 4581aa26658SKarl Rupp 4591c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 4609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xl, &xxl)); 4619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xu, &xxu)); 4622b4ed282SShri Abhyankar PetscFunctionReturn(0); 4632b4ed282SShri Abhyankar } 46492c02d66SPeter Brune 4659371c9d4SSatish Balay PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject) { 4668afaa268SBarry Smith PetscBool flg = PETSC_FALSE; 4672b4ed282SShri Abhyankar 4682b4ed282SShri Abhyankar PetscFunctionBegin; 469d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options"); 4709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL)); 4719566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL)); 4721baa6e33SBarry Smith if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL)); 473de34d3e9SBarry Smith flg = PETSC_FALSE; 4749566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_vi_monitor_residual", "Monitor residual all non-active variables; using zero for active constraints", "SNESMonitorVIResidual", flg, &flg, NULL)); 4751baa6e33SBarry Smith if (flg) PetscCall(SNESMonitorSet(snes, SNESVIMonitorResidual, PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)), NULL)); 476d0609cedSBarry Smith PetscOptionsHeadEnd(); 477ed70e96aSJungho Lee PetscFunctionReturn(0); 478ed70e96aSJungho Lee } 479