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