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