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