1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2 #include <petscdm.h> 3 4 /*@C 5 SNESVISetComputeVariableBounds - Sets a function that is called to compute the bounds on variable for 6 (differential) variable inequalities. 7 8 Input Parameters: 9 + snes - the `SNES` context 10 - compute - function that computes the bounds 11 12 Calling sequence of `compute`: 13 + snes - the `SNES` context 14 . lower - vector to hold lower bounds 15 - higher - vector to hold upper bounds 16 17 Level: advanced 18 19 Notes: 20 Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers. 21 22 For entries with no bounds you can set `PETSC_NINFINITY` or `PETSC_INFINITY` 23 24 You may use `SNESVISetVariableBounds()` to provide the bounds once if they will never change 25 26 If you have associated a `DM` with the `SNES` and provided a function to the `DM` via `DMSetVariableBounds()` that will be used automatically 27 to provide the bounds and you need not use this function. 28 29 .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `DMSetVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, 30 `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY` 31 @*/ 32 PetscErrorCode SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES snes, Vec lower, Vec higher)) 33 { 34 PetscErrorCode (*f)(SNES, PetscErrorCode (*)(SNES, Vec, Vec)); 35 36 PetscFunctionBegin; 37 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 38 PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", &f)); 39 if (f) PetscUseMethod(snes, "SNESVISetComputeVariableBounds_C", (SNES, PetscErrorCode (*)(SNES, Vec, Vec)), (snes, compute)); 40 else PetscCall(SNESVISetComputeVariableBounds_VI(snes, compute)); 41 PetscFunctionReturn(PETSC_SUCCESS); 42 } 43 44 PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES snes, SNESVIComputeVariableBoundsFn *compute) 45 { 46 PetscFunctionBegin; 47 snes->ops->computevariablebounds = compute; 48 PetscFunctionReturn(PETSC_SUCCESS); 49 } 50 51 static PetscErrorCode SNESVIMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf) 52 { 53 Vec X, F, Finactive; 54 IS isactive; 55 56 PetscFunctionBegin; 57 PetscValidHeaderSpecific(vf->viewer, PETSC_VIEWER_CLASSID, 4); 58 PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 59 PetscCall(SNESGetSolution(snes, &X)); 60 PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive)); 61 PetscCall(VecDuplicate(F, &Finactive)); 62 PetscCall(PetscObjectCompose((PetscObject)Finactive, "__Vec_bc_zero__", (PetscObject)snes)); 63 PetscCall(VecCopy(F, Finactive)); 64 PetscCall(VecISSet(Finactive, isactive, 0.0)); 65 PetscCall(ISDestroy(&isactive)); 66 PetscCall(PetscViewerPushFormat(vf->viewer, vf->format)); 67 PetscCall(VecView(Finactive, vf->viewer)); 68 PetscCall(PetscViewerPopFormat(vf->viewer)); 69 PetscCall(PetscObjectCompose((PetscObject)Finactive, "__Vec_bc_zero__", NULL)); 70 PetscCall(VecDestroy(&Finactive)); 71 PetscFunctionReturn(PETSC_SUCCESS); 72 } 73 74 static PetscErrorCode SNESVIMonitorActive(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf) 75 { 76 Vec X, F, A; 77 IS isactive; 78 79 PetscFunctionBegin; 80 PetscValidHeaderSpecific(vf->viewer, PETSC_VIEWER_CLASSID, 4); 81 PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 82 PetscCall(SNESGetSolution(snes, &X)); 83 PetscCall(SNESVIGetActiveSetIS(snes, X, F, &isactive)); 84 PetscCall(VecDuplicate(F, &A)); 85 PetscCall(PetscObjectCompose((PetscObject)A, "__Vec_bc_zero__", (PetscObject)snes)); 86 PetscCall(VecSet(A, 0.)); 87 PetscCall(VecISSet(A, isactive, 1.)); 88 PetscCall(ISDestroy(&isactive)); 89 PetscCall(PetscViewerPushFormat(vf->viewer, vf->format)); 90 PetscCall(VecView(A, vf->viewer)); 91 PetscCall(PetscViewerPopFormat(vf->viewer)); 92 PetscCall(PetscObjectCompose((PetscObject)A, "__Vec_bc_zero__", NULL)); 93 PetscCall(VecDestroy(&A)); 94 PetscFunctionReturn(PETSC_SUCCESS); 95 } 96 97 static PetscErrorCode SNESMonitorVI(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy) 98 { 99 PetscViewer viewer = (PetscViewer)dummy; 100 const PetscScalar *x, *xl, *xu, *f; 101 PetscInt i, n, act[2] = {0, 0}, fact[2], N; 102 /* Number of components that actually hit the bounds (c.f. active variables) */ 103 PetscInt act_bound[2] = {0, 0}, fact_bound[2]; 104 PetscReal rnorm, fnorm, zerotolerance = snes->vizerotolerance; 105 double tmp; 106 107 PetscFunctionBegin; 108 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4); 109 PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 110 PetscCall(VecGetSize(snes->vec_sol, &N)); 111 PetscCall(VecGetArrayRead(snes->xl, &xl)); 112 PetscCall(VecGetArrayRead(snes->xu, &xu)); 113 PetscCall(VecGetArrayRead(snes->vec_sol, &x)); 114 PetscCall(VecGetArrayRead(snes->vec_func, &f)); 115 116 rnorm = 0.0; 117 for (i = 0; i < n; i++) { 118 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]); 119 else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance && PetscRealPart(f[i]) > 0.0) act[0]++; 120 else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance && PetscRealPart(f[i]) < 0.0) act[1]++; 121 else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Can never get here"); 122 } 123 124 for (i = 0; i < n; i++) { 125 if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + zerotolerance) act_bound[0]++; 126 else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - zerotolerance) act_bound[1]++; 127 } 128 PetscCall(VecRestoreArrayRead(snes->vec_func, &f)); 129 PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 130 PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 131 PetscCall(VecRestoreArrayRead(snes->vec_sol, &x)); 132 PetscCallMPI(MPIU_Allreduce(&rnorm, &fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes))); 133 PetscCallMPI(MPIU_Allreduce(act, fact, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 134 PetscCallMPI(MPIU_Allreduce(act_bound, fact_bound, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 135 fnorm = PetscSqrtReal(fnorm); 136 137 PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel)); 138 if (snes->ntruebounds) tmp = ((double)(fact[0] + fact[1])) / ((double)snes->ntruebounds); 139 else tmp = 0.0; 140 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)); 141 142 PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel)); 143 PetscFunctionReturn(PETSC_SUCCESS); 144 } 145 146 /* 147 Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 148 || 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 149 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 150 for this trick. One assumes that the probability that W is in the null space of J is very, very small. 151 */ 152 PetscErrorCode SNESVICheckLocalMin_Private(SNES snes, Mat A, Vec F, Vec W, PetscReal fnorm, PetscBool *ismin) 153 { 154 PetscReal a1; 155 PetscBool hastranspose; 156 157 PetscFunctionBegin; 158 *ismin = PETSC_FALSE; 159 PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 160 if (hastranspose) { 161 /* Compute || J^T F|| */ 162 PetscCall(MatMultTranspose(A, F, W)); 163 PetscCall(VecNorm(W, NORM_2, &a1)); 164 PetscCall(PetscInfo(snes, "|| J^T F|| %g near zero implies found a local minimum\n", (double)(a1 / fnorm))); 165 if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE; 166 } else { 167 Vec work; 168 PetscScalar result; 169 PetscReal wnorm; 170 171 PetscCall(VecSetRandom(W, NULL)); 172 PetscCall(VecNorm(W, NORM_2, &wnorm)); 173 PetscCall(VecDuplicate(W, &work)); 174 PetscCall(MatMult(A, W, work)); 175 PetscCall(VecDot(F, work, &result)); 176 PetscCall(VecDestroy(&work)); 177 a1 = PetscAbsScalar(result) / (fnorm * wnorm); 178 PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n", (double)a1)); 179 if (a1 < 1.e-4) *ismin = PETSC_TRUE; 180 } 181 PetscFunctionReturn(PETSC_SUCCESS); 182 } 183 184 /* 185 SNESConvergedDefault_VI - Checks the convergence of the semismooth newton algorithm. 186 187 Notes: 188 The convergence criterion currently implemented is 189 merit < abstol 190 merit < rtol*merit_initial 191 */ 192 PetscErrorCode SNESConvergedDefault_VI(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gradnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 193 { 194 PetscFunctionBegin; 195 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 196 PetscAssertPointer(reason, 6); 197 198 *reason = SNES_CONVERGED_ITERATING; 199 200 if (!it) { 201 /* set parameter for default relative tolerance convergence test */ 202 snes->ttol = fnorm * snes->rtol; 203 } 204 if (fnorm != fnorm) { 205 PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n")); 206 *reason = SNES_DIVERGED_FNORM_NAN; 207 } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) { 208 PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g\n", (double)fnorm, (double)snes->abstol)); 209 *reason = SNES_CONVERGED_FNORM_ABS; 210 } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 211 PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs)); 212 *reason = SNES_DIVERGED_FUNCTION_COUNT; 213 } 214 215 if (it && !*reason) { 216 if (fnorm < snes->ttol) { 217 PetscCall(PetscInfo(snes, "Converged due to function norm %g < %g (relative tolerance)\n", (double)fnorm, (double)snes->ttol)); 218 *reason = SNES_CONVERGED_FNORM_RELATIVE; 219 } 220 } 221 PetscFunctionReturn(PETSC_SUCCESS); 222 } 223 224 /* 225 SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n. 226 227 Input Parameters: 228 . SNES - nonlinear solver context 229 230 Output Parameters: 231 . X - Bound projected X 232 233 */ 234 235 PetscErrorCode SNESVIProjectOntoBounds(SNES snes, Vec X) 236 { 237 const PetscScalar *xl, *xu; 238 PetscScalar *x; 239 PetscInt i, n; 240 241 PetscFunctionBegin; 242 PetscCall(VecGetLocalSize(X, &n)); 243 PetscCall(VecGetArray(X, &x)); 244 PetscCall(VecGetArrayRead(snes->xl, &xl)); 245 PetscCall(VecGetArrayRead(snes->xu, &xu)); 246 247 for (i = 0; i < n; i++) { 248 if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i]; 249 else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i]; 250 } 251 PetscCall(VecRestoreArray(X, &x)); 252 PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 253 PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 254 PetscFunctionReturn(PETSC_SUCCESS); 255 } 256 257 /*@ 258 SNESVIGetActiveSetIS - Gets the global indices for the active set variables 259 260 Input Parameters: 261 + snes - the `SNES` context 262 . X - the `snes` solution vector 263 - F - the nonlinear function vector 264 265 Output Parameter: 266 . ISact - active set index set 267 268 Level: developer 269 270 .seealso: [](ch_snes), `SNES`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS` 271 @*/ 272 PetscErrorCode SNESVIGetActiveSetIS(SNES snes, Vec X, Vec F, IS *ISact) 273 { 274 Vec Xl = snes->xl, Xu = snes->xu; 275 const PetscScalar *x, *f, *xl, *xu; 276 PetscInt *idx_act, i, nlocal, nloc_isact = 0, ilow, ihigh, i1 = 0; 277 PetscReal zerotolerance = snes->vizerotolerance; 278 279 PetscFunctionBegin; 280 PetscCall(VecGetLocalSize(X, &nlocal)); 281 PetscCall(VecGetOwnershipRange(X, &ilow, &ihigh)); 282 PetscCall(VecGetArrayRead(X, &x)); 283 PetscCall(VecGetArrayRead(Xl, &xl)); 284 PetscCall(VecGetArrayRead(Xu, &xu)); 285 PetscCall(VecGetArrayRead(F, &f)); 286 /* Compute active set size */ 287 for (i = 0; i < nlocal; i++) { 288 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++; 289 } 290 291 PetscCall(PetscMalloc1(nloc_isact, &idx_act)); 292 293 /* Set active set indices */ 294 for (i = 0; i < nlocal; i++) { 295 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; 296 } 297 298 /* Create active set IS */ 299 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nloc_isact, idx_act, PETSC_OWN_POINTER, ISact)); 300 301 PetscCall(VecRestoreArrayRead(X, &x)); 302 PetscCall(VecRestoreArrayRead(Xl, &xl)); 303 PetscCall(VecRestoreArrayRead(Xu, &xu)); 304 PetscCall(VecRestoreArrayRead(F, &f)); 305 PetscFunctionReturn(PETSC_SUCCESS); 306 } 307 308 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes, Vec F, Vec X, PetscReal *fnorm) 309 { 310 const PetscScalar *x, *xl, *xu, *f; 311 PetscInt i, n; 312 PetscReal rnorm, zerotolerance = snes->vizerotolerance; 313 314 PetscFunctionBegin; 315 PetscCall(VecGetLocalSize(X, &n)); 316 PetscCall(VecGetArrayRead(snes->xl, &xl)); 317 PetscCall(VecGetArrayRead(snes->xu, &xu)); 318 PetscCall(VecGetArrayRead(X, &x)); 319 PetscCall(VecGetArrayRead(F, &f)); 320 rnorm = 0.0; 321 for (i = 0; i < n; i++) { 322 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]); 323 } 324 PetscCall(VecRestoreArrayRead(F, &f)); 325 PetscCall(VecRestoreArrayRead(snes->xl, &xl)); 326 PetscCall(VecRestoreArrayRead(snes->xu, &xu)); 327 PetscCall(VecRestoreArrayRead(X, &x)); 328 PetscCallMPI(MPIU_Allreduce(&rnorm, fnorm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes))); 329 *fnorm = PetscSqrtReal(*fnorm); 330 PetscFunctionReturn(PETSC_SUCCESS); 331 } 332 333 static PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes, Vec xl, Vec xu) 334 { 335 PetscFunctionBegin; 336 PetscCall(DMComputeVariableBounds(snes->dm, xl, xu)); 337 PetscFunctionReturn(PETSC_SUCCESS); 338 } 339 340 /* 341 SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up 342 of the SNESVI nonlinear solver. 343 344 Input Parameter: 345 . snes - the SNES context 346 347 Application Interface Routine: SNESSetUp() 348 349 Notes: 350 For basic use of the SNES solvers, the user need not explicitly call 351 SNESSetUp(), since these actions will automatically occur during 352 the call to SNESSolve(). 353 */ 354 PetscErrorCode SNESSetUp_VI(SNES snes) 355 { 356 PetscInt i_start[3], i_end[3]; 357 358 PetscFunctionBegin; 359 PetscCall(SNESSetWorkVecs(snes, 1)); 360 PetscCall(SNESSetUpMatrices(snes)); 361 362 if (!snes->ops->computevariablebounds && snes->dm) { 363 PetscBool flag; 364 PetscCall(DMHasVariableBounds(snes->dm, &flag)); 365 if (flag) snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds; 366 } 367 if (!snes->usersetbounds) { 368 if (snes->ops->computevariablebounds) { 369 if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl)); 370 if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu)); 371 PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu); 372 } else if (!snes->xl && !snes->xu) { 373 /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */ 374 PetscCall(VecDuplicate(snes->vec_sol, &snes->xl)); 375 PetscCall(VecSet(snes->xl, PETSC_NINFINITY)); 376 PetscCall(VecDuplicate(snes->vec_sol, &snes->xu)); 377 PetscCall(VecSet(snes->xu, PETSC_INFINITY)); 378 } else { 379 /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 380 PetscCall(VecGetOwnershipRange(snes->vec_sol, i_start, i_end)); 381 PetscCall(VecGetOwnershipRange(snes->xl, i_start + 1, i_end + 1)); 382 PetscCall(VecGetOwnershipRange(snes->xu, i_start + 2, i_end + 2)); 383 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])) 384 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."); 385 } 386 } 387 PetscFunctionReturn(PETSC_SUCCESS); 388 } 389 PetscErrorCode SNESReset_VI(SNES snes) 390 { 391 PetscFunctionBegin; 392 PetscCall(VecDestroy(&snes->xl)); 393 PetscCall(VecDestroy(&snes->xu)); 394 snes->usersetbounds = PETSC_FALSE; 395 PetscFunctionReturn(PETSC_SUCCESS); 396 } 397 398 /* 399 SNESDestroy_VI - Destroys the private SNES_VI context that was created 400 with SNESCreate_VI(). 401 402 Input Parameter: 403 . snes - the SNES context 404 405 Application Interface Routine: SNESDestroy() 406 */ 407 PetscErrorCode SNESDestroy_VI(SNES snes) 408 { 409 PetscFunctionBegin; 410 PetscCall(PetscFree(snes->data)); 411 412 /* clear composed functions */ 413 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", NULL)); 414 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", NULL)); 415 PetscFunctionReturn(PETSC_SUCCESS); 416 } 417 418 /*@ 419 SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. `xl` <= x <= `xu`. This allows solving 420 (differential) variable inequalities. 421 422 Input Parameters: 423 + snes - the `SNES` context. 424 . xl - lower bound. 425 - xu - upper bound. 426 427 Level: advanced 428 429 Notes: 430 If this routine is not called then the lower and upper bounds are set to 431 `PETSC_NINFINITY` and `PETSC_INFINITY` respectively during `SNESSetUp()`. 432 433 Problems with bound constraints can be solved with the reduced space, `SNESVINEWTONRSLS` or semi-smooth `SNESVINEWTONSSLS` solvers. 434 435 For particular components that have no bounds you can use `PETSC_NINFINITY` or `PETSC_INFINITY` 436 437 `SNESVISetComputeVariableBounds()` can be used to provide a function that computes the bounds. This should be used if you are using, for example, grid 438 sequencing and need bounds set for a variety of vectors 439 440 .seealso: [](sec_vi), `SNES`, `SNESVIGetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY` 441 @*/ 442 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 443 { 444 PetscErrorCode (*f)(SNES, Vec, Vec); 445 446 PetscFunctionBegin; 447 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 448 PetscValidHeaderSpecific(xl, VEC_CLASSID, 2); 449 PetscValidHeaderSpecific(xu, VEC_CLASSID, 3); 450 PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESVISetVariableBounds_C", &f)); 451 if (f) PetscUseMethod(snes, "SNESVISetVariableBounds_C", (SNES, Vec, Vec), (snes, xl, xu)); 452 else PetscCall(SNESVISetVariableBounds_VI(snes, xl, xu)); 453 snes->usersetbounds = PETSC_TRUE; 454 PetscFunctionReturn(PETSC_SUCCESS); 455 } 456 457 PetscErrorCode SNESVISetVariableBounds_VI(SNES snes, Vec xl, Vec xu) 458 { 459 const PetscScalar *xxl, *xxu; 460 PetscInt i, n, cnt = 0; 461 462 PetscFunctionBegin; 463 PetscCall(SNESGetFunction(snes, &snes->vec_func, NULL, NULL)); 464 PetscCheck(snes->vec_func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first"); 465 { 466 PetscInt xlN, xuN, N; 467 PetscCall(VecGetSize(xl, &xlN)); 468 PetscCall(VecGetSize(xu, &xuN)); 469 PetscCall(VecGetSize(snes->vec_func, &N)); 470 PetscCheck(xlN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths lower bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xlN, N); 471 PetscCheck(xuN == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector lengths: upper bound = %" PetscInt_FMT " solution vector = %" PetscInt_FMT, xuN, N); 472 } 473 PetscCall(PetscObjectReference((PetscObject)xl)); 474 PetscCall(PetscObjectReference((PetscObject)xu)); 475 PetscCall(VecDestroy(&snes->xl)); 476 PetscCall(VecDestroy(&snes->xu)); 477 snes->xl = xl; 478 snes->xu = xu; 479 PetscCall(VecGetLocalSize(xl, &n)); 480 PetscCall(VecGetArrayRead(xl, &xxl)); 481 PetscCall(VecGetArrayRead(xu, &xxu)); 482 for (i = 0; i < n; i++) cnt += ((xxl[i] != PETSC_NINFINITY) || (xxu[i] != PETSC_INFINITY)); 483 484 PetscCallMPI(MPIU_Allreduce(&cnt, &snes->ntruebounds, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 485 PetscCall(VecRestoreArrayRead(xl, &xxl)); 486 PetscCall(VecRestoreArrayRead(xu, &xxu)); 487 PetscFunctionReturn(PETSC_SUCCESS); 488 } 489 490 /*@ 491 SNESVIGetVariableBounds - Gets the lower and upper bounds for the solution vector. `xl` <= x <= `xu`. These are used in solving 492 (differential) variable inequalities. 493 494 Input Parameters: 495 + snes - the `SNES` context. 496 . xl - lower bound (may be `NULL`) 497 - xu - upper bound (may be `NULL`) 498 499 Level: advanced 500 501 Note: 502 These vectors are owned by the `SNESVI` and should not be destroyed by the caller 503 504 .seealso: [](sec_vi), `SNES`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()`, `SNESVINEWTONRSLS`, `SNESVINEWTONSSLS`, `SNESSetType()`, `PETSC_NINFINITY`, `PETSC_INFINITY` 505 @*/ 506 PetscErrorCode SNESVIGetVariableBounds(SNES snes, Vec *xl, Vec *xu) 507 { 508 PetscFunctionBegin; 509 PetscCheck(snes->usersetbounds, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must set SNESVI bounds before calling SNESVIGetVariableBounds()"); 510 if (xl) *xl = snes->xl; 511 if (xu) *xu = snes->xu; 512 PetscFunctionReturn(PETSC_SUCCESS); 513 } 514 515 PetscErrorCode SNESSetFromOptions_VI(SNES snes, PetscOptionItems *PetscOptionsObject) 516 { 517 PetscBool flg = PETSC_FALSE; 518 519 PetscFunctionBegin; 520 PetscOptionsHeadBegin(PetscOptionsObject, "SNES VI options"); 521 PetscCall(PetscOptionsReal("-snes_vi_zero_tolerance", "Tolerance for considering x[] value to be on a bound", "None", snes->vizerotolerance, &snes->vizerotolerance, NULL)); 522 PetscCall(PetscOptionsBool("-snes_vi_monitor", "Monitor all non-active variables", "SNESMonitorResidual", flg, &flg, NULL)); 523 if (flg) PetscCall(SNESMonitorSet(snes, SNESMonitorVI, PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)), NULL)); 524 flg = PETSC_FALSE; 525 PetscCall(SNESMonitorSetFromOptions(snes, "-snes_vi_monitor_residual", "View residual at each iteration, using zero for active constraints", "SNESVIMonitorResidual", SNESVIMonitorResidual, NULL)); 526 PetscCall(SNESMonitorSetFromOptions(snes, "-snes_vi_monitor_active", "View active set at each iteration, using zero for inactive dofs", "SNESVIMonitorActive", SNESVIMonitorActive, NULL)); 527 PetscOptionsHeadEnd(); 528 PetscFunctionReturn(PETSC_SUCCESS); 529 } 530