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