1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 21e25c274SJed Brown #include <petscdm.h> 3d0af7cd3SBarry Smith 42176524cSBarry Smith /*@C 5f6dfbefdSBarry Smith SNESVISetComputeVariableBounds - Sets a function that is called to compute the bounds on variable for 6f6dfbefdSBarry Smith (differential) variable inequalities. 72176524cSBarry Smith 87a7aea1fSJed Brown Input parameter: 9f6dfbefdSBarry Smith + snes - the `SNES` context 10f6dfbefdSBarry Smith - compute - function that computes the bounds 11f6dfbefdSBarry Smith 12f6dfbefdSBarry Smith Calling Sequence of function: 13f6dfbefdSBarry Smith PetscErrorCode compute(SNES snes,Vec lower,Vec higher, void *ctx) 14f6dfbefdSBarry Smith 15f6dfbefdSBarry Smith + snes - the `SNES` context 16f6dfbefdSBarry Smith . lower - vector to hold lower bounds 17f6dfbefdSBarry Smith - higher - vector to hold upper bounds 182176524cSBarry Smith 192bd2b0e6SSatish Balay Level: advanced 202176524cSBarry Smith 21f6dfbefdSBarry Smith Notes: 22f0b84518SBarry Smith Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers. 23f0b84518SBarry Smith 24f6dfbefdSBarry Smith For entries with no bounds you can set `PETSC_NINFINITY` or `PETSC_INFINITY` 25c1c3a0ecSBarry Smith 26f6dfbefdSBarry Smith You may use `SNESVISetVariableBounds()` to provide the bounds once if they will never change 27f6dfbefdSBarry Smith 28f6dfbefdSBarry Smith If you have associated a `DM` with the `SNES` and provided a function to the `DM` via `DMSetVariableBounds()` that will be used automatically 29f6dfbefdSBarry Smith to provide the bounds and you need not use this function. 30f6dfbefdSBarry Smith 31*33e0caf2SBarry Smith .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `DMSetVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, 32f6dfbefdSBarry Smith 'SNESSetType()` 33aab9d709SJed Brown @*/ 34d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES, Vec, Vec)) 35d71ae5a4SJacob Faibussowitsch { 365f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(SNES, PetscErrorCode(*)(SNES, Vec, Vec)); 372176524cSBarry Smith 382176524cSBarry Smith PetscFunctionBegin; 3961589011SJed Brown PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 409566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", &f)); 41cac4c232SBarry Smith if (f) PetscUseMethod(snes, "SNESVISetComputeVariableBounds_C", (SNES, PetscErrorCode(*)(SNES, Vec, Vec)), (snes, compute)); 429566063dSJacob Faibussowitsch else PetscCall(SNESVISetComputeVariableBounds_VI(snes, compute)); 432176524cSBarry Smith PetscFunctionReturn(0); 442176524cSBarry Smith } 452176524cSBarry Smith 46d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes, SNESVIComputeVariableBoundsFunction compute) 47d71ae5a4SJacob Faibussowitsch { 4861589011SJed Brown PetscFunctionBegin; 4961589011SJed Brown snes->ops->computevariablebounds = compute; 5061589011SJed Brown PetscFunctionReturn(0); 5161589011SJed Brown } 522176524cSBarry Smith 53d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy) 54d71ae5a4SJacob Faibussowitsch { 55ffdf2a20SBarry Smith Vec X, F, Finactive; 56ffdf2a20SBarry Smith IS isactive; 57ffdf2a20SBarry Smith PetscViewer viewer = (PetscViewer)dummy; 58ffdf2a20SBarry Smith 59ffdf2a20SBarry Smith PetscFunctionBegin; 609566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 619566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &X)); 629566063dSJacob Faibussowitsch PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive)); 639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &Finactive)); 649566063dSJacob Faibussowitsch PetscCall(VecCopy(F, Finactive)); 659566063dSJacob Faibussowitsch PetscCall(VecISSet(Finactive, isactive, 0.0)); 669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isactive)); 679566063dSJacob Faibussowitsch PetscCall(VecView(Finactive, viewer)); 689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Finactive)); 69ffdf2a20SBarry Smith PetscFunctionReturn(0); 70ffdf2a20SBarry Smith } 71ffdf2a20SBarry Smith 72d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMonitorVI(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy) 73d71ae5a4SJacob Faibussowitsch { 744d4332d5SBarry Smith PetscViewer viewer = (PetscViewer)dummy; 759308c137SBarry Smith const PetscScalar *x, *xl, *xu, *f; 766fd67740SBarry Smith PetscInt i, n, act[2] = {0, 0}, fact[2], N; 776a9e2d46SJungho Lee /* Number of components that actually hit the bounds (c.f. active variables) */ 786a9e2d46SJungho Lee PetscInt act_bound[2] = {0, 0}, fact_bound[2]; 79349187a7SBarry Smith PetscReal rnorm, fnorm, zerotolerance = snes->vizerotolerance; 809d1809e2SSatish Balay double tmp; 819308c137SBarry Smith 829308c137SBarry Smith PetscFunctionBegin; 834d4332d5SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4); 849566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 859566063dSJacob Faibussowitsch PetscCall(VecGetSize(snes->vec_sol, &N)); 869566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xl, &xl)); 879566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xu, &xu)); 889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->vec_sol, &x)); 899566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->vec_func, &f)); 909308c137SBarry Smith 919308c137SBarry Smith rnorm = 0.0; 929308c137SBarry Smith for (i = 0; i < n; i++) { 93349187a7SBarry 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]); 94349187a7SBarry Smith else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++; 95349187a7SBarry Smith else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++; 96ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Can never get here"); 979308c137SBarry Smith } 986a9e2d46SJungho Lee 996a9e2d46SJungho Lee for (i = 0; i < n; i++) { 100349187a7SBarry Smith if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++; 101349187a7SBarry Smith else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++; 1026a9e2d46SJungho Lee } 1039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->vec_func, &f)); 1049566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 1059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 1069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->vec_sol, &x)); 1071c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&rnorm, &fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes))); 1081c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(act, fact, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 1091c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(act_bound, fact_bound, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 110f137e44dSBarry Smith fnorm = PetscSqrtReal(fnorm); 1116fd67740SBarry Smith 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel)); 1139d1809e2SSatish Balay if (snes->ntruebounds) tmp = ((double)(fact[0] + fact[1])) / ((double)snes->ntruebounds); 1149d1809e2SSatish Balay else tmp = 0.0; 11563a3b9bcSJacob 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)); 1166a9e2d46SJungho Lee 1179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel)); 1189308c137SBarry Smith PetscFunctionReturn(0); 1199308c137SBarry Smith } 1209308c137SBarry Smith 1212b4ed282SShri Abhyankar /* 1222b4ed282SShri Abhyankar Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 1232b4ed282SShri 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 1242b4ed282SShri Abhyankar 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 1252b4ed282SShri Abhyankar for this trick. One assumes that the probability that W is in the null space of J is very, very small. 1262b4ed282SShri Abhyankar */ 127d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVICheckLocalMin_Private(SNES snes, Mat A, Vec F, Vec W, PetscReal fnorm, PetscBool *ismin) 128d71ae5a4SJacob Faibussowitsch { 1292b4ed282SShri Abhyankar PetscReal a1; 130ace3abfcSBarry Smith PetscBool hastranspose; 1312b4ed282SShri Abhyankar 1322b4ed282SShri Abhyankar PetscFunctionBegin; 1332b4ed282SShri Abhyankar *ismin = PETSC_FALSE; 1349566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 1352b4ed282SShri Abhyankar if (hastranspose) { 1362b4ed282SShri Abhyankar /* Compute || J^T F|| */ 1379566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, F, W)); 1389566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &a1)); 1399566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "|| J^T F|| %g near zero implies found a local minimum\n", (double)(a1 / fnorm))); 1402b4ed282SShri Abhyankar if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE; 1412b4ed282SShri Abhyankar } else { 1422b4ed282SShri Abhyankar Vec work; 1432b4ed282SShri Abhyankar PetscScalar result; 1442b4ed282SShri Abhyankar PetscReal wnorm; 1452b4ed282SShri Abhyankar 1469566063dSJacob Faibussowitsch PetscCall(VecSetRandom(W, NULL)); 1479566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 1489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(W, &work)); 1499566063dSJacob Faibussowitsch PetscCall(MatMult(A, W, work)); 1509566063dSJacob Faibussowitsch PetscCall(VecDot(F, work, &result)); 1519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 1522b4ed282SShri Abhyankar a1 = PetscAbsScalar(result) / (fnorm * wnorm); 1539566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n", (double)a1)); 1542b4ed282SShri Abhyankar if (a1 < 1.e-4) *ismin = PETSC_TRUE; 1552b4ed282SShri Abhyankar } 1562b4ed282SShri Abhyankar PetscFunctionReturn(0); 1572b4ed282SShri Abhyankar } 1582b4ed282SShri Abhyankar 1592b4ed282SShri Abhyankar /* 1602b4ed282SShri Abhyankar Checks if J^T(F - J*X) = 0 1612b4ed282SShri Abhyankar */ 162d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVICheckResidual_Private(SNES snes, Mat A, Vec F, Vec X, Vec W1, Vec W2) 163d71ae5a4SJacob Faibussowitsch { 1642b4ed282SShri Abhyankar PetscReal a1, a2; 165ace3abfcSBarry Smith PetscBool hastranspose; 1662b4ed282SShri Abhyankar 1672b4ed282SShri Abhyankar PetscFunctionBegin; 1689566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 1692b4ed282SShri Abhyankar if (hastranspose) { 1709566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, W1)); 1719566063dSJacob Faibussowitsch PetscCall(VecAXPY(W1, -1.0, F)); 1722b4ed282SShri Abhyankar 1732b4ed282SShri Abhyankar /* Compute || J^T W|| */ 1749566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, W1, W2)); 1759566063dSJacob Faibussowitsch PetscCall(VecNorm(W1, NORM_2, &a1)); 1769566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &a2)); 17748a46eb9SPierre Jolivet if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n", (double)(a2 / a1))); 1782b4ed282SShri Abhyankar } 1792b4ed282SShri Abhyankar PetscFunctionReturn(0); 1802b4ed282SShri Abhyankar } 1812b4ed282SShri Abhyankar 1822b4ed282SShri Abhyankar /* 1838d359177SBarry Smith SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm. 1842b4ed282SShri Abhyankar 1852b4ed282SShri Abhyankar Notes: 1862b4ed282SShri Abhyankar The convergence criterion currently implemented is 1872b4ed282SShri Abhyankar merit < abstol 1882b4ed282SShri Abhyankar merit < rtol*merit_initial 1892b4ed282SShri Abhyankar */ 190d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESConvergedDefault_VI(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gradnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 191d71ae5a4SJacob Faibussowitsch { 1922b4ed282SShri Abhyankar PetscFunctionBegin; 1932b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1942b4ed282SShri Abhyankar PetscValidPointer(reason, 6); 1952b4ed282SShri Abhyankar 1962b4ed282SShri Abhyankar *reason = SNES_CONVERGED_ITERATING; 1972b4ed282SShri Abhyankar 1982b4ed282SShri Abhyankar if (!it) { 1992b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 2007fe79bc4SShri Abhyankar snes->ttol = fnorm * snes->rtol; 2012b4ed282SShri Abhyankar } 2027fe79bc4SShri Abhyankar if (fnorm != fnorm) { 2039566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n")); 2042b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FNORM_NAN; 2050d6f27a8SBarry Smith } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) { 2069566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g\n", (double)fnorm, (double)snes->abstol)); 2072b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_ABS; 208e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 20963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs)); 2102b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FUNCTION_COUNT; 2112b4ed282SShri Abhyankar } 2122b4ed282SShri Abhyankar 2132b4ed282SShri Abhyankar if (it && !*reason) { 2147fe79bc4SShri Abhyankar if (fnorm < snes->ttol) { 2159566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g (relative tolerance)\n", (double)fnorm, (double)snes->ttol)); 2162b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_RELATIVE; 2172b4ed282SShri Abhyankar } 2182b4ed282SShri Abhyankar } 2192b4ed282SShri Abhyankar PetscFunctionReturn(0); 2202b4ed282SShri Abhyankar } 2212b4ed282SShri Abhyankar 222c1a5e521SShri Abhyankar /* 223c1a5e521SShri Abhyankar SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n. 224c1a5e521SShri Abhyankar 225c1a5e521SShri Abhyankar Input Parameters: 226c1a5e521SShri Abhyankar . SNES - nonlinear solver context 227c1a5e521SShri Abhyankar 228c1a5e521SShri Abhyankar Output Parameters: 229c1a5e521SShri Abhyankar . X - Bound projected X 230c1a5e521SShri Abhyankar 231c1a5e521SShri Abhyankar */ 232c1a5e521SShri Abhyankar 233d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X) 234d71ae5a4SJacob Faibussowitsch { 235c1a5e521SShri Abhyankar const PetscScalar *xl, *xu; 236c1a5e521SShri Abhyankar PetscScalar *x; 237c1a5e521SShri Abhyankar PetscInt i, n; 238c1a5e521SShri Abhyankar 239c1a5e521SShri Abhyankar PetscFunctionBegin; 2409566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 2419566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 2429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xl, &xl)); 2439566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xu, &xu)); 244c1a5e521SShri Abhyankar 245c1a5e521SShri Abhyankar for (i = 0; i < n; i++) { 24610f5ae6bSBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i]; 24710f5ae6bSBarry Smith else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i]; 248c1a5e521SShri Abhyankar } 2499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 2509566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 2519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 252c1a5e521SShri Abhyankar PetscFunctionReturn(0); 253c1a5e521SShri Abhyankar } 254c1a5e521SShri Abhyankar 255693ddf92SShri Abhyankar /* 256693ddf92SShri Abhyankar SNESVIGetActiveSetIndices - Gets the global indices for the active set variables 257693ddf92SShri Abhyankar 2587a7aea1fSJed Brown Input parameter: 259693ddf92SShri Abhyankar . snes - the SNES context 260693ddf92SShri Abhyankar . X - the snes solution vector 261693ddf92SShri Abhyankar . F - the nonlinear function vector 262693ddf92SShri Abhyankar 2637a7aea1fSJed Brown Output parameter: 264693ddf92SShri Abhyankar . ISact - active set index set 265693ddf92SShri Abhyankar */ 266d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact) 267d71ae5a4SJacob Faibussowitsch { 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; 2749566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nlocal)); 2759566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh)); 2769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xl, &xl)); 2789566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xu, &xu)); 2799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &f)); 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 2859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc_isact, &idx_act)); 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 */ 2939566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact)); 294d950fb63SShri Abhyankar 2959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xl, &xl)); 2979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xu, &xu)); 2989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &f)); 299d950fb63SShri Abhyankar PetscFunctionReturn(0); 300d950fb63SShri Abhyankar } 301d950fb63SShri Abhyankar 302d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVICreateIndexSets_RS(SNES snes, Vec X, Vec F, IS *ISact, IS *ISinact) 303d71ae5a4SJacob Faibussowitsch { 304077aedafSJed Brown PetscInt rstart, rend; 305693ddf92SShri Abhyankar 306693ddf92SShri Abhyankar PetscFunctionBegin; 3079566063dSJacob Faibussowitsch PetscCall(SNESVIGetActiveSetIS(snes, X, F, ISact)); 3089566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &rstart, &rend)); 3099566063dSJacob Faibussowitsch PetscCall(ISComplement(*ISact, rstart, rend, ISinact)); 310693ddf92SShri Abhyankar PetscFunctionReturn(0); 311693ddf92SShri Abhyankar } 312693ddf92SShri Abhyankar 313d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm) 314d71ae5a4SJacob Faibussowitsch { 315fe835674SShri Abhyankar const PetscScalar *x, *xl, *xu, *f; 316fe835674SShri Abhyankar PetscInt i, n; 317feb237baSPierre Jolivet PetscReal rnorm, zerotolerance = snes->vizerotolerance; 318fe835674SShri Abhyankar 319fe835674SShri Abhyankar PetscFunctionBegin; 3209566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 3219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xl, &xl)); 3229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xu, &xu)); 3239566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 3249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &f)); 325fe835674SShri Abhyankar rnorm = 0.0; 326fe835674SShri Abhyankar for (i = 0; i < n; i++) { 327349187a7SBarry 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]); 3288f918527SKarl Rupp } 3299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &f)); 3309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 3319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 3329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 3331c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes))); 33462d1f40fSBarry Smith *fnorm = PetscSqrtReal(*fnorm); 335fe835674SShri Abhyankar PetscFunctionReturn(0); 336fe835674SShri Abhyankar } 337fe835674SShri Abhyankar 338d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu) 339d71ae5a4SJacob Faibussowitsch { 34008da532bSDmitry Karpeev PetscFunctionBegin; 3419566063dSJacob Faibussowitsch PetscCall(DMComputeVariableBounds(snes->dm, xl, xu)); 34208da532bSDmitry Karpeev PetscFunctionReturn(0); 34308da532bSDmitry Karpeev } 34408da532bSDmitry Karpeev 3452b4ed282SShri Abhyankar /* 346c2fc9fa9SBarry Smith SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up 3472b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 3482b4ed282SShri Abhyankar 3492b4ed282SShri Abhyankar Input Parameter: 3502b4ed282SShri Abhyankar . snes - the SNES context 3512b4ed282SShri Abhyankar 3522b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 3532b4ed282SShri Abhyankar 3542b4ed282SShri Abhyankar Notes: 3552b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 3562b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 3572b4ed282SShri Abhyankar the call to SNESSolve(). 3582b4ed282SShri Abhyankar */ 359d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_VI(SNES snes) 360d71ae5a4SJacob Faibussowitsch { 3612b4ed282SShri Abhyankar PetscInt i_start[3], i_end[3]; 3622b4ed282SShri Abhyankar 3632b4ed282SShri Abhyankar PetscFunctionBegin; 3649566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 1)); 3659566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 3662b4ed282SShri Abhyankar 36708da532bSDmitry Karpeev if (!snes->ops->computevariablebounds && snes->dm) { 368a201590fSDmitry Karpeev PetscBool flag; 3699566063dSJacob Faibussowitsch PetscCall(DMHasVariableBounds(snes->dm, &flag)); 370ad540459SPierre Jolivet if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds; 37136b2a9f3SBarry Smith } 372a201590fSDmitry Karpeev if (!snes->usersetbounds) { 373c2fc9fa9SBarry Smith if (snes->ops->computevariablebounds) { 3749566063dSJacob Faibussowitsch if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl)); 3759566063dSJacob Faibussowitsch if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu)); 376dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu); 3771aa26658SKarl Rupp } else if (!snes->xl && !snes->xu) { 3782176524cSBarry Smith /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */ 3799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &snes->xl)); 3809566063dSJacob Faibussowitsch PetscCall(VecSet(snes->xl, PETSC_NINFINITY)); 3819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &snes->xu)); 3829566063dSJacob Faibussowitsch PetscCall(VecSet(snes->xu, PETSC_INFINITY)); 3832b4ed282SShri Abhyankar } else { 3842b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 3859566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end)); 3869566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1)); 3879566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2)); 3882b4ed282SShri 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])) 3892b4ed282SShri 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."); 3902b4ed282SShri Abhyankar } 391a201590fSDmitry Karpeev } 3922b4ed282SShri Abhyankar PetscFunctionReturn(0); 3932b4ed282SShri Abhyankar } 394d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_VI(SNES snes) 395d71ae5a4SJacob Faibussowitsch { 3962176524cSBarry Smith PetscFunctionBegin; 3979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xl)); 3989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xu)); 3992d6615e8SDmitry Karpeev snes->usersetbounds = PETSC_FALSE; 4002176524cSBarry Smith PetscFunctionReturn(0); 4012176524cSBarry Smith } 4022176524cSBarry Smith 4032b4ed282SShri Abhyankar /* 4042b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 4052b4ed282SShri Abhyankar with SNESCreate_VI(). 4062b4ed282SShri Abhyankar 4072b4ed282SShri Abhyankar Input Parameter: 4082b4ed282SShri Abhyankar . snes - the SNES context 4092b4ed282SShri Abhyankar 4102b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 4112b4ed282SShri Abhyankar */ 412d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_VI(SNES snes) 413d71ae5a4SJacob Faibussowitsch { 4142b4ed282SShri Abhyankar PetscFunctionBegin; 4159566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 4162b4ed282SShri Abhyankar 4172b4ed282SShri Abhyankar /* clear composed functions */ 4182e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL)); 4192e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL)); 4202b4ed282SShri Abhyankar PetscFunctionReturn(0); 4212b4ed282SShri Abhyankar } 4227fe79bc4SShri Abhyankar 4235559b345SBarry Smith /*@ 424f0b84518SBarry Smith SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving 425f6dfbefdSBarry Smith (differential) variable inequalities. 4262b4ed282SShri Abhyankar 4272b4ed282SShri Abhyankar Input Parameters: 428f6dfbefdSBarry Smith + snes - the `SNES` context. 4292b4ed282SShri Abhyankar . xl - lower bound. 430a2b725a8SWilliam Gropp - xu - upper bound. 4312b4ed282SShri Abhyankar 4322b4ed282SShri Abhyankar Notes: 4332b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 434f6dfbefdSBarry Smith `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`. 43584c105d7SBarry Smith 436f0b84518SBarry Smith Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers. 437f0b84518SBarry Smith 438f6dfbefdSBarry Smith For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY` 439f6dfbefdSBarry Smith 440*33e0caf2SBarry Smith `SNESVISetComputeVariableBounds()` can be used to provide a function that computes the bounds. This should be used if you are using, for example, grid 441f6dfbefdSBarry Smith sequencing and need bounds set for a variety of vectors 442f6dfbefdSBarry Smith 4432bd2b0e6SSatish Balay Level: advanced 4442bd2b0e6SSatish Balay 445*33e0caf2SBarry Smith .seealso: [](sec_vi), `SNES`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, SNESVINEWTONRSLS, SNESVINEWTONSSLS, 'SNESSetType()` 4465559b345SBarry Smith @*/ 447d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 448d71ae5a4SJacob Faibussowitsch { 4495f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(SNES, Vec, Vec); 4502b4ed282SShri Abhyankar 4512b4ed282SShri Abhyankar PetscFunctionBegin; 4522b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 4532b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl, VEC_CLASSID, 2); 4542b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu, VEC_CLASSID, 3); 4559566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f)); 456cac4c232SBarry Smith if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu)); 4579566063dSJacob Faibussowitsch else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu)); 458a201590fSDmitry Karpeev snes->usersetbounds = PETSC_TRUE; 45961589011SJed Brown PetscFunctionReturn(0); 46061589011SJed Brown } 46161589011SJed Brown 462d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu) 463d71ae5a4SJacob Faibussowitsch { 46461589011SJed Brown const PetscScalar *xxl, *xxu; 46561589011SJed Brown PetscInt i, n, cnt = 0; 46661589011SJed Brown 46761589011SJed Brown PetscFunctionBegin; 4689566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL)); 4695f80ce2aSJacob Faibussowitsch PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first"); 470077aedafSJed Brown { 471077aedafSJed Brown PetscInt xlN, xuN, N; 4729566063dSJacob Faibussowitsch PetscCall(VecGetSize(xl, &xlN)); 4739566063dSJacob Faibussowitsch PetscCall(VecGetSize(xu, &xuN)); 4749566063dSJacob Faibussowitsch PetscCall(VecGetSize(snes->vec_func, &N)); 47563a3b9bcSJacob Faibussowitsch PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N); 47663a3b9bcSJacob Faibussowitsch PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N); 477077aedafSJed Brown } 4789566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xl)); 4799566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xu)); 4809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xl)); 4819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xu)); 482c2fc9fa9SBarry Smith snes->xl = xl; 483c2fc9fa9SBarry Smith snes->xu = xu; 4849566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xl, &n)); 4859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xl, &xxl)); 4869566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xu, &xxu)); 487e270355aSBarry Smith for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY)); 4881aa26658SKarl Rupp 4891c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 4909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xl, &xxl)); 4919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xu, &xxu)); 4922b4ed282SShri Abhyankar PetscFunctionReturn(0); 4932b4ed282SShri Abhyankar } 49492c02d66SPeter Brune 495d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject) 496d71ae5a4SJacob Faibussowitsch { 4978afaa268SBarry Smith PetscBool flg = PETSC_FALSE; 4982b4ed282SShri Abhyankar 4992b4ed282SShri Abhyankar PetscFunctionBegin; 500d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options"); 5019566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL)); 5029566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL)); 5031baa6e33SBarry Smith if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL)); 504de34d3e9SBarry Smith flg = PETSC_FALSE; 5059566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_vi_monitor_residual", "Monitor residual all non-active variables; using zero for active constraints", "SNESMonitorVIResidual", flg, &flg, NULL)); 5061baa6e33SBarry Smith if (flg) PetscCall(SNESMonitorSet(snes, SNESVIMonitorResidual, PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)), NULL)); 507d0609cedSBarry Smith PetscOptionsHeadEnd(); 508ed70e96aSJungho Lee PetscFunctionReturn(0); 509ed70e96aSJungho Lee } 510