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