1 #include "viimpl.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 100 PetscFunctionBegin; 101 ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr); 102 ierr = VecGetSize(snes->vec_sol,&N);CHKERRQ(ierr); 103 ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 104 ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 105 ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 106 ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr); 107 108 rnorm = 0.0; 109 for (i=0; i<n; i++) { 110 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]); 111 else if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8 && PetscRealPart(f[i]) >= 0.0) act[0]++; 112 else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8 && PetscRealPart(f[i]) <= 0.0) act[1]++; 113 else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Can never get here"); 114 } 115 116 for (i=0; i<n; i++) { 117 if (PetscRealPart(x[i]) <= PetscRealPart(xl[i]) + 1.e-8) act_bound[0]++; 118 else if (PetscRealPart(x[i]) >= PetscRealPart(xu[i]) - 1.e-8) act_bound[1]++; 119 } 120 ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr); 121 ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 122 ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 123 ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 124 ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr); 125 ierr = MPI_Allreduce(act,fact,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr); 126 ierr = MPI_Allreduce(act_bound,fact_bound,2,MPIU_INT,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr); 127 fnorm = PetscSqrtReal(fnorm); 128 129 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 130 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),((double)(fact[0]+fact[1]))/((double)snes->ntruebounds));CHKERRQ(ierr); 131 132 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 133 PetscFunctionReturn(0); 134 } 135 136 /* 137 Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 138 || 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 139 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 140 for this trick. One assumes that the probability that W is in the null space of J is very, very small. 141 */ 142 #undef __FUNCT__ 143 #define __FUNCT__ "SNESVICheckLocalMin_Private" 144 PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 145 { 146 PetscReal a1; 147 PetscErrorCode ierr; 148 PetscBool hastranspose; 149 150 PetscFunctionBegin; 151 *ismin = PETSC_FALSE; 152 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 153 if (hastranspose) { 154 /* Compute || J^T F|| */ 155 ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 156 ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 157 ierr = PetscInfo1(snes,"|| J^T F|| %g near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr); 158 if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 159 } else { 160 Vec work; 161 PetscScalar result; 162 PetscReal wnorm; 163 164 ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 165 ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 166 ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 167 ierr = MatMult(A,W,work);CHKERRQ(ierr); 168 ierr = VecDot(F,work,&result);CHKERRQ(ierr); 169 ierr = VecDestroy(&work);CHKERRQ(ierr); 170 a1 = PetscAbsScalar(result)/(fnorm*wnorm); 171 ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr); 172 if (a1 < 1.e-4) *ismin = PETSC_TRUE; 173 } 174 PetscFunctionReturn(0); 175 } 176 177 /* 178 Checks if J^T(F - J*X) = 0 179 */ 180 #undef __FUNCT__ 181 #define __FUNCT__ "SNESVICheckResidual_Private" 182 PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 183 { 184 PetscReal a1,a2; 185 PetscErrorCode ierr; 186 PetscBool hastranspose; 187 188 PetscFunctionBegin; 189 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 190 if (hastranspose) { 191 ierr = MatMult(A,X,W1);CHKERRQ(ierr); 192 ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 193 194 /* Compute || J^T W|| */ 195 ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 196 ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 197 ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 198 if (a1 != 0.0) { 199 ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr); 200 } 201 } 202 PetscFunctionReturn(0); 203 } 204 205 /* 206 SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm. 207 208 Notes: 209 The convergence criterion currently implemented is 210 merit < abstol 211 merit < rtol*merit_initial 212 */ 213 #undef __FUNCT__ 214 #define __FUNCT__ "SNESDefaultConverged_VI" 215 PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 216 { 217 PetscErrorCode ierr; 218 219 PetscFunctionBegin; 220 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 221 PetscValidPointer(reason,6); 222 223 *reason = SNES_CONVERGED_ITERATING; 224 225 if (!it) { 226 /* set parameter for default relative tolerance convergence test */ 227 snes->ttol = fnorm*snes->rtol; 228 } 229 if (fnorm != fnorm) { 230 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 231 *reason = SNES_DIVERGED_FNORM_NAN; 232 } else if (fnorm < snes->abstol) { 233 ierr = PetscInfo2(snes,"Converged due to function norm %g < %g\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr); 234 *reason = SNES_CONVERGED_FNORM_ABS; 235 } else if (snes->nfuncs >= snes->max_funcs) { 236 ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 237 *reason = SNES_DIVERGED_FUNCTION_COUNT; 238 } 239 240 if (it && !*reason) { 241 if (fnorm < snes->ttol) { 242 ierr = PetscInfo2(snes,"Converged due to function norm %g < %g (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr); 243 *reason = SNES_CONVERGED_FNORM_RELATIVE; 244 } 245 } 246 PetscFunctionReturn(0); 247 } 248 249 250 /* -------------------------------------------------------------------------- */ 251 /* 252 SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n. 253 254 Input Parameters: 255 . SNES - nonlinear solver context 256 257 Output Parameters: 258 . X - Bound projected X 259 260 */ 261 262 #undef __FUNCT__ 263 #define __FUNCT__ "SNESVIProjectOntoBounds" 264 PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X) 265 { 266 PetscErrorCode ierr; 267 const PetscScalar *xl,*xu; 268 PetscScalar *x; 269 PetscInt i,n; 270 271 PetscFunctionBegin; 272 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 273 ierr = VecGetArray(X,&x);CHKERRQ(ierr); 274 ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 275 ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 276 277 for(i = 0;i<n;i++) { 278 if (PetscRealPart(x[i]) < PetscRealPart(xl[i])) x[i] = xl[i]; 279 else if (PetscRealPart(x[i]) > PetscRealPart(xu[i])) x[i] = xu[i]; 280 } 281 ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 282 ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 283 ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 284 PetscFunctionReturn(0); 285 } 286 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "SNESVIGetActiveSetIS" 290 /* 291 SNESVIGetActiveSetIndices - Gets the global indices for the active set variables 292 293 Input parameter 294 . snes - the SNES context 295 . X - the snes solution vector 296 . F - the nonlinear function vector 297 298 Output parameter 299 . ISact - active set index set 300 */ 301 PetscErrorCode SNESVIGetActiveSetIS(SNES snes,Vec X,Vec F,IS* ISact) 302 { 303 PetscErrorCode ierr; 304 Vec Xl=snes->xl,Xu=snes->xu; 305 const PetscScalar *x,*f,*xl,*xu; 306 PetscInt *idx_act,i,nlocal,nloc_isact=0,ilow,ihigh,i1=0; 307 308 PetscFunctionBegin; 309 ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 310 ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr); 311 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 312 ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr); 313 ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr); 314 ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 315 /* Compute active set size */ 316 for (i=0; i < nlocal;i++) { 317 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++; 318 } 319 320 ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 321 322 /* Set active set indices */ 323 for(i=0; i < nlocal; i++) { 324 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; 325 } 326 327 /* Create active set IS */ 328 ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_OWN_POINTER,ISact);CHKERRQ(ierr); 329 330 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 331 ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr); 332 ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr); 333 ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 334 PetscFunctionReturn(0); 335 } 336 337 #undef __FUNCT__ 338 #define __FUNCT__ "SNESVICreateIndexSets_RS" 339 PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec F,IS* ISact,IS* ISinact) 340 { 341 PetscErrorCode ierr; 342 343 PetscFunctionBegin; 344 ierr = SNESVIGetActiveSetIS(snes,X,F,ISact);CHKERRQ(ierr); 345 ierr = ISComplement(*ISact,X->map->rstart,X->map->rend,ISinact);CHKERRQ(ierr); 346 PetscFunctionReturn(0); 347 } 348 349 #undef __FUNCT__ 350 #define __FUNCT__ "SNESVIComputeInactiveSetFnorm" 351 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X, PetscReal *fnorm) 352 { 353 PetscErrorCode ierr; 354 const PetscScalar *x,*xl,*xu,*f; 355 PetscInt i,n; 356 PetscReal rnorm; 357 358 PetscFunctionBegin; 359 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 360 ierr = VecGetArrayRead(snes->xl,&xl);CHKERRQ(ierr); 361 ierr = VecGetArrayRead(snes->xu,&xu);CHKERRQ(ierr); 362 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 363 ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 364 rnorm = 0.0; 365 for (i=0; i<n; i++) { 366 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]); 367 } 368 ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 369 ierr = VecRestoreArrayRead(snes->xl,&xl);CHKERRQ(ierr); 370 ierr = VecRestoreArrayRead(snes->xu,&xu);CHKERRQ(ierr); 371 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 372 ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr); 373 *fnorm = PetscSqrtReal(*fnorm); 374 PetscFunctionReturn(0); 375 } 376 377 #undef __FUNCT__ 378 #define __FUNCT__ "SNESVIDMComputeVariableBounds" 379 PetscErrorCode SNESVIDMComputeVariableBounds(SNES snes,Vec xl, Vec xu) 380 { 381 PetscErrorCode ierr; 382 PetscFunctionBegin; 383 ierr = DMComputeVariableBounds(snes->dm, xl, xu); CHKERRQ(ierr); 384 PetscFunctionReturn(0); 385 } 386 387 388 /* -------------------------------------------------------------------------- */ 389 /* 390 SNESSetUp_VI - Does setup common to all VI solvers -- basically makes sure bounds have been properly set up 391 of the SNESVI nonlinear solver. 392 393 Input Parameter: 394 . snes - the SNES context 395 . x - the solution vector 396 397 Application Interface Routine: SNESSetUp() 398 399 Notes: 400 For basic use of the SNES solvers, the user need not explicitly call 401 SNESSetUp(), since these actions will automatically occur during 402 the call to SNESSolve(). 403 */ 404 #undef __FUNCT__ 405 #define __FUNCT__ "SNESSetUp_VI" 406 PetscErrorCode SNESSetUp_VI(SNES snes) 407 { 408 PetscErrorCode ierr; 409 PetscInt i_start[3],i_end[3]; 410 411 PetscFunctionBegin; 412 413 ierr = SNESDefaultGetWork(snes,1);CHKERRQ(ierr); 414 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 415 416 if(!snes->ops->computevariablebounds && snes->dm) { 417 snes->ops->computevariablebounds = SNESVIDMComputeVariableBounds; 418 } 419 if (snes->ops->computevariablebounds) { 420 if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);} 421 if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);} 422 ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr); 423 } else if (!snes->xl && !snes->xu) { 424 /* If the lower and upper bound on variables are not set, set it to -Inf and Inf */ 425 ierr = VecDuplicate(snes->vec_sol, &snes->xl);CHKERRQ(ierr); 426 ierr = VecSet(snes->xl,SNES_VI_NINF);CHKERRQ(ierr); 427 ierr = VecDuplicate(snes->vec_sol, &snes->xu);CHKERRQ(ierr); 428 ierr = VecSet(snes->xu,SNES_VI_INF);CHKERRQ(ierr); 429 } else { 430 /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 431 ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 432 ierr = VecGetOwnershipRange(snes->xl,i_start+1,i_end+1);CHKERRQ(ierr); 433 ierr = VecGetOwnershipRange(snes->xu,i_start+2,i_end+2);CHKERRQ(ierr); 434 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])) 435 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."); 436 } 437 438 PetscFunctionReturn(0); 439 } 440 /* -------------------------------------------------------------------------- */ 441 #undef __FUNCT__ 442 #define __FUNCT__ "SNESReset_VI" 443 PetscErrorCode SNESReset_VI(SNES snes) 444 { 445 PetscErrorCode ierr; 446 447 PetscFunctionBegin; 448 ierr = VecDestroy(&snes->xl);CHKERRQ(ierr); 449 ierr = VecDestroy(&snes->xu);CHKERRQ(ierr); 450 PetscFunctionReturn(0); 451 } 452 453 /* 454 SNESDestroy_VI - Destroys the private SNES_VI context that was created 455 with SNESCreate_VI(). 456 457 Input Parameter: 458 . snes - the SNES context 459 460 Application Interface Routine: SNESDestroy() 461 */ 462 #undef __FUNCT__ 463 #define __FUNCT__ "SNESDestroy_VI" 464 PetscErrorCode SNESDestroy_VI(SNES snes) 465 { 466 PetscErrorCode ierr; 467 468 PetscFunctionBegin; 469 ierr = PetscFree(snes->data);CHKERRQ(ierr); 470 471 /* clear composed functions */ 472 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 473 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 474 PetscFunctionReturn(0); 475 } 476 477 /* -------------------------------------------------------------------------- */ 478 479 /* 480 481 These line searches are common for all the VI solvers 482 */ 483 extern PetscErrorCode SNESSolve_VISS(SNES); 484 485 #undef __FUNCT__ 486 #define __FUNCT__ "SNESLineSearchNo_VI" 487 488 /* 489 This routine does not actually do a line search but it takes a full newton 490 step while ensuring that the new iterates remain within the constraints. 491 492 */ 493 PetscErrorCode SNESLineSearchNo_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 494 { 495 PetscErrorCode ierr; 496 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 497 498 PetscFunctionBegin; 499 *flag = PETSC_TRUE; 500 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 501 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 502 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 503 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 504 if (snes->ops->postcheckstep) { 505 ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 506 } 507 if (changed_y) { 508 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 509 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 510 } 511 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 512 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 513 if (!snes->domainerror) { 514 if (snes->ops->solve != SNESSolve_VISS) { 515 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 516 } else { 517 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 518 } 519 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 520 } 521 if (snes->ls_monitor) { 522 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 523 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr); 524 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 525 } 526 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 527 PetscFunctionReturn(0); 528 } 529 530 /* -------------------------------------------------------------------------- */ 531 #undef __FUNCT__ 532 #define __FUNCT__ "SNESLineSearchNoNorms_VI" 533 534 /* 535 This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 536 */ 537 PetscErrorCode SNESLineSearchNoNorms_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 538 { 539 PetscErrorCode ierr; 540 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 541 542 PetscFunctionBegin; 543 *flag = PETSC_TRUE; 544 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 545 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 546 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 547 if (snes->ops->postcheckstep) { 548 ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 549 } 550 if (changed_y) { 551 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 552 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 553 } 554 555 /* don't evaluate function the last time through */ 556 if (snes->iter < snes->max_its-1) { 557 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 558 } 559 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 560 PetscFunctionReturn(0); 561 } 562 563 /* -------------------------------------------------------------------------- */ 564 565 #undef __FUNCT__ 566 #define __FUNCT__ "SNESLineSearchCubic_VI" 567 /* 568 This routine implements a cubic line search while doing a projection on the variable bounds 569 */ 570 PetscErrorCode SNESLineSearchCubic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 571 { 572 PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 573 PetscReal minlambda,lambda,lambdatemp; 574 #if defined(PETSC_USE_COMPLEX) 575 PetscScalar cinitslope; 576 #endif 577 PetscErrorCode ierr; 578 PetscInt count; 579 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 580 MPI_Comm comm; 581 582 PetscFunctionBegin; 583 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 584 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 585 *flag = PETSC_TRUE; 586 587 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 588 if (*ynorm == 0.0) { 589 if (snes->ls_monitor) { 590 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 591 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 592 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 593 } 594 *gnorm = fnorm; 595 ierr = VecCopy(x,w);CHKERRQ(ierr); 596 ierr = VecCopy(f,g);CHKERRQ(ierr); 597 *flag = PETSC_FALSE; 598 goto theend1; 599 } 600 if (*ynorm > snes->maxstep) { /* Step too big, so scale back */ 601 if (snes->ls_monitor) { 602 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 603 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: Scaling step by %g old ynorm %g\n",(double)(snes->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr); 604 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 605 } 606 ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr); 607 *ynorm = snes->maxstep; 608 } 609 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 610 minlambda = snes->steptol/rellength; 611 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 612 #if defined(PETSC_USE_COMPLEX) 613 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 614 initslope = PetscRealPart(cinitslope); 615 #else 616 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 617 #endif 618 if (initslope > 0.0) initslope = -initslope; 619 if (initslope == 0.0) initslope = -1.0; 620 621 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 622 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 623 if (snes->nfuncs >= snes->max_funcs) { 624 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 625 *flag = PETSC_FALSE; 626 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 627 goto theend1; 628 } 629 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 630 if (snes->ops->solve != SNESSolve_VISS) { 631 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 632 } else { 633 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 634 } 635 if (snes->domainerror) { 636 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 637 PetscFunctionReturn(0); 638 } 639 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 640 ierr = PetscInfo4(snes,"Initial fnorm %g gnorm %g alpha %g initslope %g\n",(double)fnorm,(double)*gnorm,(double)snes->ls_alpha,(double)initslope);CHKERRQ(ierr); 641 if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* Sufficient reduction */ 642 if (snes->ls_monitor) { 643 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 644 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr); 645 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 646 } 647 goto theend1; 648 } 649 650 /* Fit points with quadratic */ 651 lambda = 1.0; 652 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 653 lambdaprev = lambda; 654 gnormprev = *gnorm; 655 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 656 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 657 else lambda = lambdatemp; 658 659 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 660 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 661 if (snes->nfuncs >= snes->max_funcs) { 662 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 663 *flag = PETSC_FALSE; 664 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 665 goto theend1; 666 } 667 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 668 if (snes->ops->solve != SNESSolve_VISS) { 669 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 670 } else { 671 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 672 } 673 if (snes->domainerror) { 674 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 675 PetscFunctionReturn(0); 676 } 677 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 678 if (snes->ls_monitor) { 679 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 680 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: gnorm after quadratic fit %g\n",(double)(*gnorm));CHKERRQ(ierr); 681 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 682 } 683 if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm ) { /* sufficient reduction */ 684 if (snes->ls_monitor) { 685 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 686 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 687 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 688 } 689 goto theend1; 690 } 691 692 /* Fit points with cubic */ 693 count = 1; 694 while (PETSC_TRUE) { 695 if (lambda <= minlambda) { 696 if (snes->ls_monitor) { 697 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 698 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 699 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)minlambda,(double)lambda,(double)initslope);CHKERRQ(ierr); 700 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 701 } 702 *flag = PETSC_FALSE; 703 break; 704 } 705 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 706 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 707 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 708 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 709 d = b*b - 3*a*initslope; 710 if (d < 0.0) d = 0.0; 711 if (a == 0.0) { 712 lambdatemp = -initslope/(2.0*b); 713 } else { 714 lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 715 } 716 lambdaprev = lambda; 717 gnormprev = *gnorm; 718 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 719 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 720 else lambda = lambdatemp; 721 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 722 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 723 if (snes->nfuncs >= snes->max_funcs) { 724 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 725 ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)lambda,(double)initslope);CHKERRQ(ierr); 726 *flag = PETSC_FALSE; 727 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 728 break; 729 } 730 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 731 if (snes->ops->solve != SNESSolve_VISS) { 732 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 733 } else { 734 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 735 } 736 if (snes->domainerror) { 737 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 738 PetscFunctionReturn(0); 739 } 740 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 741 if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* is reduction enough? */ 742 if (snes->ls_monitor) { 743 ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %g lambda=%18.16e\n",(double)(*gnorm),(double)lambda);CHKERRQ(ierr); 744 } 745 break; 746 } else { 747 if (snes->ls_monitor) { 748 ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %g lambda=%18.16e\n",(double)(*gnorm),(double)lambda);CHKERRQ(ierr); 749 } 750 } 751 count++; 752 } 753 theend1: 754 /* Optional user-defined check for line search step validity */ 755 if (snes->ops->postcheckstep && *flag) { 756 ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 757 if (changed_y) { 758 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 759 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 760 } 761 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 762 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 763 if (snes->ops->solve != SNESSolve_VISS) { 764 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 765 } else { 766 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 767 } 768 if (snes->domainerror) { 769 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 770 PetscFunctionReturn(0); 771 } 772 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 773 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 774 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 775 } 776 } 777 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 778 PetscFunctionReturn(0); 779 } 780 781 #undef __FUNCT__ 782 #define __FUNCT__ "SNESLineSearchQuadratic_VI" 783 /* 784 This routine does a quadratic line search while keeping the iterates within the variable bounds 785 */ 786 PetscErrorCode SNESLineSearchQuadratic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 787 { 788 /* 789 Note that for line search purposes we work with with the related 790 minimization problem: 791 min z(x): R^n -> R, 792 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 793 */ 794 PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 795 #if defined(PETSC_USE_COMPLEX) 796 PetscScalar cinitslope; 797 #endif 798 PetscErrorCode ierr; 799 PetscInt count; 800 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 801 802 PetscFunctionBegin; 803 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 804 *flag = PETSC_TRUE; 805 806 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 807 if (*ynorm == 0.0) { 808 if (snes->ls_monitor) { 809 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 810 ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 811 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 812 } 813 *gnorm = fnorm; 814 ierr = VecCopy(x,w);CHKERRQ(ierr); 815 ierr = VecCopy(f,g);CHKERRQ(ierr); 816 *flag = PETSC_FALSE; 817 goto theend2; 818 } 819 if (*ynorm > snes->maxstep) { /* Step too big, so scale back */ 820 ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr); 821 *ynorm = snes->maxstep; 822 } 823 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 824 minlambda = snes->steptol/rellength; 825 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 826 #if defined(PETSC_USE_COMPLEX) 827 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 828 initslope = PetscRealPart(cinitslope); 829 #else 830 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 831 #endif 832 if (initslope > 0.0) initslope = -initslope; 833 if (initslope == 0.0) initslope = -1.0; 834 ierr = PetscInfo1(snes,"Initslope %g \n",(double)initslope);CHKERRQ(ierr); 835 836 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 837 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 838 if (snes->nfuncs >= snes->max_funcs) { 839 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 840 *flag = PETSC_FALSE; 841 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 842 goto theend2; 843 } 844 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 845 if (snes->ops->solve != SNESSolve_VISS) { 846 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 847 } else { 848 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 849 } 850 if (snes->domainerror) { 851 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 852 PetscFunctionReturn(0); 853 } 854 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 855 if ((*gnorm)*(*gnorm) <= (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* Sufficient reduction */ 856 if (snes->ls_monitor) { 857 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 858 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: Using full step: fnorm %g gnorm %g\n",(double)fnorm,(double)(*gnorm));CHKERRQ(ierr); 859 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 860 } 861 goto theend2; 862 } 863 864 /* Fit points with quadratic */ 865 lambda = 1.0; 866 count = 1; 867 while (PETSC_TRUE) { 868 if (lambda <= minlambda) { /* bad luck; use full step */ 869 if (snes->ls_monitor) { 870 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 871 ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 872 ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)lambda,(double)initslope);CHKERRQ(ierr); 873 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 874 } 875 ierr = VecCopy(x,w);CHKERRQ(ierr); 876 *flag = PETSC_FALSE; 877 break; 878 } 879 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 880 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 881 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 882 else lambda = lambdatemp; 883 884 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 885 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 886 if (snes->nfuncs >= snes->max_funcs) { 887 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 888 ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)(*gnorm),(double)(*ynorm),(double)lambda,(double)initslope);CHKERRQ(ierr); 889 *flag = PETSC_FALSE; 890 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 891 break; 892 } 893 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 894 if (snes->domainerror) { 895 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 896 PetscFunctionReturn(0); 897 } 898 if (snes->ops->solve != SNESSolve_VISS) { 899 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 900 } else { 901 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 902 } 903 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 904 if ((*gnorm)*(*gnorm) < (1.0 - snes->ls_alpha)*fnorm*fnorm) { /* sufficient reduction */ 905 if (snes->ls_monitor) { 906 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 907 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line Search: Quadratically determined step, lambda=%g\n",(double)lambda);CHKERRQ(ierr); 908 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 909 } 910 break; 911 } 912 count++; 913 } 914 theend2: 915 /* Optional user-defined check for line search step validity */ 916 if (snes->ops->postcheckstep) { 917 ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 918 if (changed_y) { 919 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 920 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 921 } 922 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 923 ierr = SNESComputeFunction(snes,w,g); 924 if (snes->domainerror) { 925 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 926 PetscFunctionReturn(0); 927 } 928 if (snes->ops->solve != SNESSolve_VISS) { 929 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 930 } else { 931 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 932 } 933 934 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 935 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 936 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 937 } 938 } 939 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 940 PetscFunctionReturn(0); 941 } 942 943 944 #undef __FUNCT__ 945 #define __FUNCT__ "SNESVISetVariableBounds" 946 /*@ 947 SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 948 949 Input Parameters: 950 . snes - the SNES context. 951 . xl - lower bound. 952 . xu - upper bound. 953 954 Notes: 955 If this routine is not called then the lower and upper bounds are set to 956 SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp(). 957 958 Level: advanced 959 960 @*/ 961 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 962 { 963 PetscErrorCode ierr,(*f)(SNES,Vec,Vec); 964 965 PetscFunctionBegin; 966 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 967 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 968 PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 969 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESVISetVariableBounds_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr); 970 if (!f) {ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr);} 971 ierr = PetscUseMethod(snes,"SNESVISetVariableBounds_C",(SNES,Vec,Vec),(snes,xl,xu));CHKERRQ(ierr); 972 PetscFunctionReturn(0); 973 } 974 975 EXTERN_C_BEGIN 976 #undef __FUNCT__ 977 #define __FUNCT__ "SNESVISetVariableBounds_VI" 978 PetscErrorCode SNESVISetVariableBounds_VI(SNES snes,Vec xl,Vec xu) 979 { 980 PetscErrorCode ierr; 981 const PetscScalar *xxl,*xxu; 982 PetscInt i,n, cnt = 0; 983 984 PetscFunctionBegin; 985 ierr = SNESGetFunction(snes,&snes->vec_func,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 986 if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first"); 987 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); 988 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); 989 ierr = SNESSetType(snes,SNESVIRS);CHKERRQ(ierr); 990 ierr = PetscObjectReference((PetscObject)xl);CHKERRQ(ierr); 991 ierr = PetscObjectReference((PetscObject)xu);CHKERRQ(ierr); 992 ierr = VecDestroy(&snes->xl);CHKERRQ(ierr); 993 ierr = VecDestroy(&snes->xu);CHKERRQ(ierr); 994 snes->xl = xl; 995 snes->xu = xu; 996 ierr = VecGetLocalSize(xl,&n);CHKERRQ(ierr); 997 ierr = VecGetArrayRead(xl,&xxl);CHKERRQ(ierr); 998 ierr = VecGetArrayRead(xu,&xxu);CHKERRQ(ierr); 999 for (i=0; i<n; i++) { 1000 cnt += ((xxl[i] != SNES_VI_NINF) || (xxu[i] != SNES_VI_INF)); 1001 } 1002 ierr = MPI_Allreduce(&cnt,&snes->ntruebounds,1,MPIU_INT,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr); 1003 ierr = VecRestoreArrayRead(xl,&xxl);CHKERRQ(ierr); 1004 ierr = VecRestoreArrayRead(xu,&xxu);CHKERRQ(ierr); 1005 PetscFunctionReturn(0); 1006 } 1007 EXTERN_C_END 1008 1009 EXTERN_C_BEGIN 1010 #undef __FUNCT__ 1011 #define __FUNCT__ "SNESLineSearchSetType_VI" 1012 PetscErrorCode SNESLineSearchSetType_VI(SNES snes, SNESLineSearchType type) 1013 { 1014 PetscErrorCode ierr; 1015 PetscFunctionBegin; 1016 1017 switch (type) { 1018 case SNES_LS_BASIC: 1019 ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 1020 break; 1021 case SNES_LS_BASIC_NONORMS: 1022 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 1023 break; 1024 case SNES_LS_QUADRATIC: 1025 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 1026 break; 1027 case SNES_LS_CUBIC: 1028 ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 1029 break; 1030 default: 1031 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 1032 break; 1033 } 1034 snes->ls_type = type; 1035 PetscFunctionReturn(0); 1036 } 1037 EXTERN_C_END 1038 1039 #undef __FUNCT__ 1040 #define __FUNCT__ "SNESSetFromOptions_VI" 1041 PetscErrorCode SNESSetFromOptions_VI(SNES snes) 1042 { 1043 PetscErrorCode ierr; 1044 PetscBool flg; 1045 1046 PetscFunctionBegin; 1047 ierr = PetscOptionsHead("SNES VI options");CHKERRQ(ierr); 1048 ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 1049 if (flg) { 1050 ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr); 1051 } 1052 ierr = PetscOptionsTail();CHKERRQ(ierr); 1053 PetscFunctionReturn(0); 1054 } 1055