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