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 1220f4b53cSBarry Smith Calling Sequence of `compute`: 1320f4b53cSBarry Smith $ PetscErrorCode compute(SNES snes, Vec lower, Vec higher) 14f6dfbefdSBarry Smith + snes - the `SNES` context 15f6dfbefdSBarry Smith . lower - vector to hold lower bounds 16f6dfbefdSBarry Smith - higher - vector to hold upper bounds 172176524cSBarry Smith 182bd2b0e6SSatish Balay Level: advanced 192176524cSBarry Smith 20f6dfbefdSBarry Smith Notes: 21f0b84518SBarry Smith Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers. 22f0b84518SBarry Smith 23f6dfbefdSBarry Smith For entries with no bounds you can set `PETSC_NINFINITY` or `PETSC_INFINITY` 24c1c3a0ecSBarry Smith 25f6dfbefdSBarry Smith You may use `SNESVISetVariableBounds()` to provide the bounds once if they will never change 26f6dfbefdSBarry Smith 27f6dfbefdSBarry Smith If you have associated a `DM` with the `SNES` and provided a function to the `DM` via `DMSetVariableBounds()` that will be used automatically 28f6dfbefdSBarry Smith to provide the bounds and you need not use this function. 29f6dfbefdSBarry Smith 3033e0caf2SBarry Smith .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `DMSetVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, 31f6dfbefdSBarry Smith 'SNESSetType()` 32aab9d709SJed Brown @*/ 33d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES, Vec, Vec)) 34d71ae5a4SJacob Faibussowitsch { 355f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(SNES, PetscErrorCode (*)(SNES, Vec, Vec)); 362176524cSBarry Smith 372176524cSBarry Smith PetscFunctionBegin; 3861589011SJed Brown PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 399566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", &f)); 40cac4c232SBarry Smith if (f) PetscUseMethod(snes, "SNESVISetComputeVariableBounds_C", (SNES, PetscErrorCode(*)(SNES, Vec, Vec)), (snes, compute)); 419566063dSJacob Faibussowitsch else PetscCall(SNESVISetComputeVariableBounds_VI(snes, compute)); 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 432176524cSBarry Smith } 442176524cSBarry Smith 45d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes, SNESVIComputeVariableBoundsFunction compute) 46d71ae5a4SJacob Faibussowitsch { 4761589011SJed Brown PetscFunctionBegin; 4861589011SJed Brown snes->ops->computevariablebounds = compute; 493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5061589011SJed Brown } 512176524cSBarry Smith 52d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy) 53d71ae5a4SJacob Faibussowitsch { 54ffdf2a20SBarry Smith Vec X, F, Finactive; 55ffdf2a20SBarry Smith IS isactive; 56ffdf2a20SBarry Smith PetscViewer viewer = (PetscViewer)dummy; 57ffdf2a20SBarry Smith 58ffdf2a20SBarry Smith PetscFunctionBegin; 599566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 609566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &X)); 619566063dSJacob Faibussowitsch PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive)); 629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &Finactive)); 639566063dSJacob Faibussowitsch PetscCall(VecCopy(F, Finactive)); 649566063dSJacob Faibussowitsch PetscCall(VecISSet(Finactive, isactive, 0.0)); 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isactive)); 669566063dSJacob Faibussowitsch PetscCall(VecView(Finactive, viewer)); 679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Finactive)); 683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69ffdf2a20SBarry Smith } 70ffdf2a20SBarry Smith 71d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMonitorVI(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy) 72d71ae5a4SJacob Faibussowitsch { 734d4332d5SBarry Smith PetscViewer viewer = (PetscViewer)dummy; 749308c137SBarry Smith const PetscScalar *x, *xl, *xu, *f; 756fd67740SBarry Smith PetscInt i, n, act[2] = {0, 0}, fact[2], N; 766a9e2d46SJungho Lee /* Number of components that actually hit the bounds (c.f. active variables) */ 776a9e2d46SJungho Lee PetscInt act_bound[2] = {0, 0}, fact_bound[2]; 78349187a7SBarry Smith PetscReal rnorm, fnorm, zerotolerance = snes->vizerotolerance; 799d1809e2SSatish Balay double tmp; 809308c137SBarry Smith 819308c137SBarry Smith PetscFunctionBegin; 824d4332d5SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4); 839566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 849566063dSJacob Faibussowitsch PetscCall(VecGetSize(snes->vec_sol, &N)); 859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xl, &xl)); 869566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xu, &xu)); 879566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->vec_sol, &x)); 889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->vec_func, &f)); 899308c137SBarry Smith 909308c137SBarry Smith rnorm = 0.0; 919308c137SBarry Smith for (i = 0; i < n; i++) { 92349187a7SBarry 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]); 93349187a7SBarry Smith else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++; 94349187a7SBarry Smith else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++; 95ce94432eSBarry Smith else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Can never get here"); 969308c137SBarry Smith } 976a9e2d46SJungho Lee 986a9e2d46SJungho Lee for (i = 0; i < n; i++) { 99349187a7SBarry Smith if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++; 100349187a7SBarry Smith else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++; 1016a9e2d46SJungho Lee } 1029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->vec_func, &f)); 1039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 1049566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 1059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->vec_sol, &x)); 1061c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&rnorm, &fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes))); 1071c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(act, fact, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 1081c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(act_bound, fact_bound, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 109f137e44dSBarry Smith fnorm = PetscSqrtReal(fnorm); 1106fd67740SBarry Smith 1119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel)); 1129d1809e2SSatish Balay if (snes->ntruebounds) tmp = ((double)(fact[0] + fact[1])) / ((double)snes->ntruebounds); 1139d1809e2SSatish Balay else tmp = 0.0; 11463a3b9bcSJacob 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)); 1156a9e2d46SJungho Lee 1169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel)); 1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1189308c137SBarry Smith } 1199308c137SBarry Smith 1202b4ed282SShri Abhyankar /* 1212b4ed282SShri Abhyankar Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 1222b4ed282SShri 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 1232b4ed282SShri Abhyankar 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 1242b4ed282SShri Abhyankar for this trick. One assumes that the probability that W is in the null space of J is very, very small. 1252b4ed282SShri Abhyankar */ 126d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVICheckLocalMin_Private(SNES snes, Mat A, Vec F, Vec W, PetscReal fnorm, PetscBool *ismin) 127d71ae5a4SJacob Faibussowitsch { 1282b4ed282SShri Abhyankar PetscReal a1; 129ace3abfcSBarry Smith PetscBool hastranspose; 1302b4ed282SShri Abhyankar 1312b4ed282SShri Abhyankar PetscFunctionBegin; 1322b4ed282SShri Abhyankar *ismin = PETSC_FALSE; 1339566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 1342b4ed282SShri Abhyankar if (hastranspose) { 1352b4ed282SShri Abhyankar /* Compute || J^T F|| */ 1369566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, F, W)); 1379566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &a1)); 1389566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "|| J^T F|| %g near zero implies found a local minimum\n", (double)(a1 / fnorm))); 1392b4ed282SShri Abhyankar if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE; 1402b4ed282SShri Abhyankar } else { 1412b4ed282SShri Abhyankar Vec work; 1422b4ed282SShri Abhyankar PetscScalar result; 1432b4ed282SShri Abhyankar PetscReal wnorm; 1442b4ed282SShri Abhyankar 1459566063dSJacob Faibussowitsch PetscCall(VecSetRandom(W, NULL)); 1469566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 1479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(W, &work)); 1489566063dSJacob Faibussowitsch PetscCall(MatMult(A, W, work)); 1499566063dSJacob Faibussowitsch PetscCall(VecDot(F, work, &result)); 1509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 1512b4ed282SShri Abhyankar a1 = PetscAbsScalar(result) / (fnorm * wnorm); 1529566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n", (double)a1)); 1532b4ed282SShri Abhyankar if (a1 < 1.e-4) *ismin = PETSC_TRUE; 1542b4ed282SShri Abhyankar } 1553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1562b4ed282SShri Abhyankar } 1572b4ed282SShri Abhyankar 1582b4ed282SShri Abhyankar /* 1592b4ed282SShri Abhyankar Checks if J^T(F - J*X) = 0 1602b4ed282SShri Abhyankar */ 161d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVICheckResidual_Private(SNES snes, Mat A, Vec F, Vec X, Vec W1, Vec W2) 162d71ae5a4SJacob Faibussowitsch { 1632b4ed282SShri Abhyankar PetscReal a1, a2; 164ace3abfcSBarry Smith PetscBool hastranspose; 1652b4ed282SShri Abhyankar 1662b4ed282SShri Abhyankar PetscFunctionBegin; 1679566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 1682b4ed282SShri Abhyankar if (hastranspose) { 1699566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, W1)); 1709566063dSJacob Faibussowitsch PetscCall(VecAXPY(W1, -1.0, F)); 1712b4ed282SShri Abhyankar 1722b4ed282SShri Abhyankar /* Compute || J^T W|| */ 1739566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, W1, W2)); 1749566063dSJacob Faibussowitsch PetscCall(VecNorm(W1, NORM_2, &a1)); 1759566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &a2)); 17648a46eb9SPierre Jolivet if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n", (double)(a2 / a1))); 1772b4ed282SShri Abhyankar } 1783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1792b4ed282SShri Abhyankar } 1802b4ed282SShri Abhyankar 1812b4ed282SShri Abhyankar /* 1828d359177SBarry Smith SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm. 1832b4ed282SShri Abhyankar 1842b4ed282SShri Abhyankar Notes: 1852b4ed282SShri Abhyankar The convergence criterion currently implemented is 1862b4ed282SShri Abhyankar merit < abstol 1872b4ed282SShri Abhyankar merit < rtol*merit_initial 1882b4ed282SShri Abhyankar */ 189d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESConvergedDefault_VI(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gradnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 190d71ae5a4SJacob Faibussowitsch { 1912b4ed282SShri Abhyankar PetscFunctionBegin; 1922b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1932b4ed282SShri Abhyankar PetscValidPointer(reason, 6); 1942b4ed282SShri Abhyankar 1952b4ed282SShri Abhyankar *reason = SNES_CONVERGED_ITERATING; 1962b4ed282SShri Abhyankar 1972b4ed282SShri Abhyankar if (!it) { 1982b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 1997fe79bc4SShri Abhyankar snes->ttol = fnorm * snes->rtol; 2002b4ed282SShri Abhyankar } 2017fe79bc4SShri Abhyankar if (fnorm != fnorm) { 2029566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n")); 2032b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FNORM_NAN; 2040d6f27a8SBarry Smith } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) { 2059566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g\n", (double)fnorm, (double)snes->abstol)); 2062b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_ABS; 207e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 20863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs)); 2092b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FUNCTION_COUNT; 2102b4ed282SShri Abhyankar } 2112b4ed282SShri Abhyankar 2122b4ed282SShri Abhyankar if (it && !*reason) { 2137fe79bc4SShri Abhyankar if (fnorm < snes->ttol) { 2149566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g (relative tolerance)\n", (double)fnorm, (double)snes->ttol)); 2152b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_RELATIVE; 2162b4ed282SShri Abhyankar } 2172b4ed282SShri Abhyankar } 2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2192b4ed282SShri Abhyankar } 2202b4ed282SShri Abhyankar 221c1a5e521SShri Abhyankar /* 222c1a5e521SShri Abhyankar SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n. 223c1a5e521SShri Abhyankar 224c1a5e521SShri Abhyankar Input Parameters: 225c1a5e521SShri Abhyankar . SNES - nonlinear solver context 226c1a5e521SShri Abhyankar 227c1a5e521SShri Abhyankar Output Parameters: 228c1a5e521SShri Abhyankar . X - Bound projected X 229c1a5e521SShri Abhyankar 230c1a5e521SShri Abhyankar */ 231c1a5e521SShri Abhyankar 232d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X) 233d71ae5a4SJacob Faibussowitsch { 234c1a5e521SShri Abhyankar const PetscScalar *xl, *xu; 235c1a5e521SShri Abhyankar PetscScalar *x; 236c1a5e521SShri Abhyankar PetscInt i, n; 237c1a5e521SShri Abhyankar 238c1a5e521SShri Abhyankar PetscFunctionBegin; 2399566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 2409566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 2419566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xl, &xl)); 2429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xu, &xu)); 243c1a5e521SShri Abhyankar 244c1a5e521SShri Abhyankar for (i = 0; i < n; i++) { 24510f5ae6bSBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i]; 24610f5ae6bSBarry Smith else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i]; 247c1a5e521SShri Abhyankar } 2489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 2499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 2509566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 252c1a5e521SShri Abhyankar } 253c1a5e521SShri Abhyankar 254693ddf92SShri Abhyankar /* 255693ddf92SShri Abhyankar SNESVIGetActiveSetIndices - Gets the global indices for the active set variables 256693ddf92SShri Abhyankar 2577a7aea1fSJed Brown Input parameter: 258693ddf92SShri Abhyankar . snes - the SNES context 259693ddf92SShri Abhyankar . X - the snes solution vector 260693ddf92SShri Abhyankar . F - the nonlinear function vector 261693ddf92SShri Abhyankar 2627a7aea1fSJed Brown Output parameter: 263693ddf92SShri Abhyankar . ISact - active set index set 264693ddf92SShri Abhyankar */ 265d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact) 266d71ae5a4SJacob Faibussowitsch { 267c2fc9fa9SBarry Smith Vec Xl = snes->xl, Xu = snes->xu; 268693ddf92SShri Abhyankar const PetscScalar *x, *f, *xl, *xu; 269693ddf92SShri Abhyankar PetscInt *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0; 270349187a7SBarry Smith PetscReal zerotolerance = snes->vizerotolerance; 271d950fb63SShri Abhyankar 272d950fb63SShri Abhyankar PetscFunctionBegin; 2739566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nlocal)); 2749566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh)); 2759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xl, &xl)); 2779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xu, &xu)); 2789566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &f)); 279693ddf92SShri Abhyankar /* Compute active set size */ 280d950fb63SShri Abhyankar for (i = 0; i < nlocal; i++) { 281349187a7SBarry 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++; 282d950fb63SShri Abhyankar } 283693ddf92SShri Abhyankar 2849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc_isact, &idx_act)); 285d950fb63SShri Abhyankar 286693ddf92SShri Abhyankar /* Set active set indices */ 287d950fb63SShri Abhyankar for (i = 0; i < nlocal; i++) { 288349187a7SBarry 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; 289d950fb63SShri Abhyankar } 290d950fb63SShri Abhyankar 291693ddf92SShri Abhyankar /* Create active set IS */ 2929566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact)); 293d950fb63SShri Abhyankar 2949566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xl, &xl)); 2969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xu, &xu)); 2979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &f)); 2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 299d950fb63SShri Abhyankar } 300d950fb63SShri Abhyankar 301d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVICreateIndexSets_RS(SNES snes, Vec X, Vec F, IS *ISact, IS *ISinact) 302d71ae5a4SJacob Faibussowitsch { 303077aedafSJed Brown PetscInt rstart, rend; 304693ddf92SShri Abhyankar 305693ddf92SShri Abhyankar PetscFunctionBegin; 3069566063dSJacob Faibussowitsch PetscCall(SNESVIGetActiveSetIS(snes, X, F, ISact)); 3079566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &rstart, &rend)); 3089566063dSJacob Faibussowitsch PetscCall(ISComplement(*ISact, rstart, rend, ISinact)); 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 310693ddf92SShri Abhyankar } 311693ddf92SShri Abhyankar 312d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm) 313d71ae5a4SJacob Faibussowitsch { 314fe835674SShri Abhyankar const PetscScalar *x, *xl, *xu, *f; 315fe835674SShri Abhyankar PetscInt i, n; 316feb237baSPierre Jolivet PetscReal rnorm, zerotolerance = snes->vizerotolerance; 317fe835674SShri Abhyankar 318fe835674SShri Abhyankar PetscFunctionBegin; 3199566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 3209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xl, &xl)); 3219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(snes->xu, &xu)); 3229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 3239566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &f)); 324fe835674SShri Abhyankar rnorm = 0.0; 325fe835674SShri Abhyankar for (i = 0; i < n; i++) { 326349187a7SBarry 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]); 3278f918527SKarl Rupp } 3289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &f)); 3299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 3309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 3319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 3321c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes))); 33362d1f40fSBarry Smith *fnorm = PetscSqrtReal(*fnorm); 3343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 335fe835674SShri Abhyankar } 336fe835674SShri Abhyankar 337d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu) 338d71ae5a4SJacob Faibussowitsch { 33908da532bSDmitry Karpeev PetscFunctionBegin; 3409566063dSJacob Faibussowitsch PetscCall(DMComputeVariableBounds(snes->dm, xl, xu)); 3413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34208da532bSDmitry Karpeev } 34308da532bSDmitry Karpeev 3442b4ed282SShri Abhyankar /* 345c2fc9fa9SBarry Smith SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up 3462b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 3472b4ed282SShri Abhyankar 3482b4ed282SShri Abhyankar Input Parameter: 3492b4ed282SShri Abhyankar . snes - the SNES context 3502b4ed282SShri Abhyankar 3512b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 3522b4ed282SShri Abhyankar 3532b4ed282SShri Abhyankar Notes: 3542b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 3552b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 3562b4ed282SShri Abhyankar the call to SNESSolve(). 3572b4ed282SShri Abhyankar */ 358d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_VI(SNES snes) 359d71ae5a4SJacob Faibussowitsch { 3602b4ed282SShri Abhyankar PetscInt i_start[3], i_end[3]; 3612b4ed282SShri Abhyankar 3622b4ed282SShri Abhyankar PetscFunctionBegin; 3639566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 1)); 3649566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 3652b4ed282SShri Abhyankar 36608da532bSDmitry Karpeev if (!snes->ops->computevariablebounds && snes->dm) { 367a201590fSDmitry Karpeev PetscBool flag; 3689566063dSJacob Faibussowitsch PetscCall(DMHasVariableBounds(snes->dm, &flag)); 369ad540459SPierre Jolivet if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds; 37036b2a9f3SBarry Smith } 371a201590fSDmitry Karpeev if (!snes->usersetbounds) { 372c2fc9fa9SBarry Smith if (snes->ops->computevariablebounds) { 3739566063dSJacob Faibussowitsch if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl)); 3749566063dSJacob Faibussowitsch if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu)); 375dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu); 3761aa26658SKarl Rupp } else if (!snes->xl && !snes->xu) { 3772176524cSBarry Smith /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */ 3789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &snes->xl)); 3799566063dSJacob Faibussowitsch PetscCall(VecSet(snes->xl, PETSC_NINFINITY)); 3809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &snes->xu)); 3819566063dSJacob Faibussowitsch PetscCall(VecSet(snes->xu, PETSC_INFINITY)); 3822b4ed282SShri Abhyankar } else { 3832b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 3849566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end)); 3859566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1)); 3869566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2)); 3872b4ed282SShri 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])) 3882b4ed282SShri 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."); 3892b4ed282SShri Abhyankar } 390a201590fSDmitry Karpeev } 3913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3922b4ed282SShri Abhyankar } 393d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_VI(SNES snes) 394d71ae5a4SJacob Faibussowitsch { 3952176524cSBarry Smith PetscFunctionBegin; 3969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xl)); 3979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xu)); 3982d6615e8SDmitry Karpeev snes->usersetbounds = PETSC_FALSE; 3993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4002176524cSBarry Smith } 4012176524cSBarry Smith 4022b4ed282SShri Abhyankar /* 4032b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 4042b4ed282SShri Abhyankar with SNESCreate_VI(). 4052b4ed282SShri Abhyankar 4062b4ed282SShri Abhyankar Input Parameter: 4072b4ed282SShri Abhyankar . snes - the SNES context 4082b4ed282SShri Abhyankar 4092b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 4102b4ed282SShri Abhyankar */ 411d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_VI(SNES snes) 412d71ae5a4SJacob Faibussowitsch { 4132b4ed282SShri Abhyankar PetscFunctionBegin; 4149566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 4152b4ed282SShri Abhyankar 4162b4ed282SShri Abhyankar /* clear composed functions */ 4172e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL)); 4182e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL)); 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4202b4ed282SShri Abhyankar } 4217fe79bc4SShri Abhyankar 4225559b345SBarry Smith /*@ 423f0b84518SBarry Smith SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving 424f6dfbefdSBarry Smith (differential) variable inequalities. 4252b4ed282SShri Abhyankar 4262b4ed282SShri Abhyankar Input Parameters: 427f6dfbefdSBarry Smith + snes - the `SNES` context. 4282b4ed282SShri Abhyankar . xl - lower bound. 429a2b725a8SWilliam Gropp - xu - upper bound. 4302b4ed282SShri Abhyankar 431*2fe279fdSBarry Smith Level: advanced 432*2fe279fdSBarry Smith 4332b4ed282SShri Abhyankar Notes: 4342b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 435f6dfbefdSBarry Smith `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`. 43684c105d7SBarry Smith 437f0b84518SBarry Smith Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers. 438f0b84518SBarry Smith 439f6dfbefdSBarry Smith For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY` 440f6dfbefdSBarry Smith 44133e0caf2SBarry Smith `SNESVISetComputeVariableBounds()` can be used to provide a function that computes the bounds. This should be used if you are using, for example, grid 442f6dfbefdSBarry Smith sequencing and need bounds set for a variety of vectors 443f6dfbefdSBarry Smith 444cf836535SBlaise Bourdin .seealso: [](sec_vi), `SNES`, `SNESVIGetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, 'SNESSetType()` 4455559b345SBarry Smith @*/ 446d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 447d71ae5a4SJacob Faibussowitsch { 4485f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(SNES, Vec, Vec); 4492b4ed282SShri Abhyankar 4502b4ed282SShri Abhyankar PetscFunctionBegin; 4512b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 4522b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl, VEC_CLASSID, 2); 4532b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu, VEC_CLASSID, 3); 4549566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f)); 455cac4c232SBarry Smith if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu)); 4569566063dSJacob Faibussowitsch else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu)); 457a201590fSDmitry Karpeev snes->usersetbounds = PETSC_TRUE; 4583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45961589011SJed Brown } 46061589011SJed Brown 461d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu) 462d71ae5a4SJacob Faibussowitsch { 46361589011SJed Brown const PetscScalar *xxl, *xxu; 46461589011SJed Brown PetscInt i, n, cnt = 0; 46561589011SJed Brown 46661589011SJed Brown PetscFunctionBegin; 4679566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL)); 4685f80ce2aSJacob Faibussowitsch PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first"); 469077aedafSJed Brown { 470077aedafSJed Brown PetscInt xlN, xuN, N; 4719566063dSJacob Faibussowitsch PetscCall(VecGetSize(xl, &xlN)); 4729566063dSJacob Faibussowitsch PetscCall(VecGetSize(xu, &xuN)); 4739566063dSJacob Faibussowitsch PetscCall(VecGetSize(snes->vec_func, &N)); 47463a3b9bcSJacob Faibussowitsch PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N); 47563a3b9bcSJacob Faibussowitsch PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N); 476077aedafSJed Brown } 4779566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xl)); 4789566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)xu)); 4799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xl)); 4809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&snes->xu)); 481c2fc9fa9SBarry Smith snes->xl = xl; 482c2fc9fa9SBarry Smith snes->xu = xu; 4839566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xl, &n)); 4849566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xl, &xxl)); 4859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xu, &xxu)); 486e270355aSBarry Smith for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY)); 4871aa26658SKarl Rupp 4881c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 4899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xl, &xxl)); 4909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xu, &xxu)); 4913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4922b4ed282SShri Abhyankar } 49392c02d66SPeter Brune 494cf836535SBlaise Bourdin /*@ 495cf836535SBlaise Bourdin SNESVIGetVariableBounds - Gets the lower and upper bounds for the solution vector. xl <= x <= xu. This allows solving 496cf836535SBlaise Bourdin (differential) variable inequalities. 497cf836535SBlaise Bourdin 498cf836535SBlaise Bourdin Input Parameters: 499cf836535SBlaise Bourdin + snes - the `SNES` context. 500cf836535SBlaise Bourdin . xl - lower bound (may be `NULL`) 501cf836535SBlaise Bourdin - xu - upper bound (may be `NULL`) 502cf836535SBlaise Bourdin 503*2fe279fdSBarry Smith Level: advanced 504*2fe279fdSBarry Smith 505cf836535SBlaise Bourdin Notes: 506cf836535SBlaise Bourdin These vectors are owned by the `SNESVI` and should not be destroyed by the caller 507cf836535SBlaise Bourdin 508cf836535SBlaise Bourdin .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, SNESVINEWTONRSLS, SNESVINEWTONSSLS, 'SNESSetType()` 509cf836535SBlaise Bourdin @*/ 510cf836535SBlaise Bourdin PetscErrorCode SNESVIGetVariableBounds(SNES snes, Vec *xl, Vec *xu) 511cf836535SBlaise Bourdin { 512cf836535SBlaise Bourdin PetscFunctionBegin; 513cf836535SBlaise Bourdin PetscCheck(snes->usersetbounds, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must set SNESVI bounds before calling SNESVIGetVariableBounds()"); 514cf836535SBlaise Bourdin if (xl) *xl = snes->xl; 515cf836535SBlaise Bourdin if (xu) *xu = snes->xu; 5163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 517cf836535SBlaise Bourdin } 518cf836535SBlaise Bourdin 519d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject) 520d71ae5a4SJacob Faibussowitsch { 5218afaa268SBarry Smith PetscBool flg = PETSC_FALSE; 5222b4ed282SShri Abhyankar 5232b4ed282SShri Abhyankar PetscFunctionBegin; 524d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options"); 5259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL)); 5269566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL)); 5271baa6e33SBarry Smith if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL)); 528de34d3e9SBarry Smith flg = PETSC_FALSE; 5299566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_vi_monitor_residual", "Monitor residual all non-active variables; using zero for active constraints", "SNESMonitorVIResidual", flg, &flg, NULL)); 5301baa6e33SBarry Smith if (flg) PetscCall(SNESMonitorSet(snes, SNESVIMonitorResidual, PETSC_VIEWER_DRAW_(PetscObjectComm((PetscObject)snes)), NULL)); 531d0609cedSBarry Smith PetscOptionsHeadEnd(); 5323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 533ed70e96aSJungho Lee } 534