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