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