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