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