1 #define PETSCSNES_DLL 2 3 #include "../src/snes/impls/vi/viimpl.h" /*I "petscsnes.h" I*/ 4 #include "../include/private/kspimpl.h" 5 6 #undef __FUNCT__ 7 #define __FUNCT__ "SNESMonitorVI" 8 PetscErrorCode SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 9 { 10 PetscErrorCode ierr; 11 SNES_VI *vi = (SNES_VI*)snes->data; 12 PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy; 13 const PetscScalar *x,*xl,*xu,*f; 14 PetscInt i,n,act = 0; 15 PetscReal rnorm,fnorm; 16 17 PetscFunctionBegin; 18 if (!dummy) { 19 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 20 } 21 ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr); 22 ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr); 23 ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr); 24 ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 25 ierr = VecGetArrayRead(snes->vec_func,&f);CHKERRQ(ierr); 26 27 rnorm = 0.0; 28 for (i=0; i<n; i++) { 29 if (((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) rnorm += f[i]*f[i]; 30 else act++; 31 } 32 ierr = VecRestoreArrayRead(snes->vec_func,&f);CHKERRQ(ierr); 33 ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr); 34 ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr); 35 ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 36 ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr); 37 fnorm = sqrt(fnorm); 38 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES VI Function norm %14.12e Active constraints %D\n",its,fnorm,act);CHKERRQ(ierr); 39 if (!dummy) { 40 ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 41 } 42 PetscFunctionReturn(0); 43 } 44 45 /* 46 Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 47 || 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 48 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 49 for this trick. One assumes that the probability that W is in the null space of J is very, very small. 50 */ 51 #undef __FUNCT__ 52 #define __FUNCT__ "SNESVICheckLocalMin_Private" 53 PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 54 { 55 PetscReal a1; 56 PetscErrorCode ierr; 57 PetscBool hastranspose; 58 59 PetscFunctionBegin; 60 *ismin = PETSC_FALSE; 61 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 62 if (hastranspose) { 63 /* Compute || J^T F|| */ 64 ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 65 ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 66 ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr); 67 if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 68 } else { 69 Vec work; 70 PetscScalar result; 71 PetscReal wnorm; 72 73 ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 74 ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 75 ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 76 ierr = MatMult(A,W,work);CHKERRQ(ierr); 77 ierr = VecDot(F,work,&result);CHKERRQ(ierr); 78 ierr = VecDestroy(work);CHKERRQ(ierr); 79 a1 = PetscAbsScalar(result)/(fnorm*wnorm); 80 ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr); 81 if (a1 < 1.e-4) *ismin = PETSC_TRUE; 82 } 83 PetscFunctionReturn(0); 84 } 85 86 /* 87 Checks if J^T(F - J*X) = 0 88 */ 89 #undef __FUNCT__ 90 #define __FUNCT__ "SNESVICheckResidual_Private" 91 PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 92 { 93 PetscReal a1,a2; 94 PetscErrorCode ierr; 95 PetscBool hastranspose; 96 97 PetscFunctionBegin; 98 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 99 if (hastranspose) { 100 ierr = MatMult(A,X,W1);CHKERRQ(ierr); 101 ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 102 103 /* Compute || J^T W|| */ 104 ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 105 ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 106 ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 107 if (a1 != 0.0) { 108 ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 109 } 110 } 111 PetscFunctionReturn(0); 112 } 113 114 /* 115 SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm. 116 117 Notes: 118 The convergence criterion currently implemented is 119 merit < abstol 120 merit < rtol*merit_initial 121 */ 122 #undef __FUNCT__ 123 #define __FUNCT__ "SNESDefaultConverged_VI" 124 PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 125 { 126 PetscErrorCode ierr; 127 128 PetscFunctionBegin; 129 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 130 PetscValidPointer(reason,6); 131 132 *reason = SNES_CONVERGED_ITERATING; 133 134 if (!it) { 135 /* set parameter for default relative tolerance convergence test */ 136 snes->ttol = fnorm*snes->rtol; 137 } 138 if (fnorm != fnorm) { 139 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 140 *reason = SNES_DIVERGED_FNORM_NAN; 141 } else if (fnorm < snes->abstol) { 142 ierr = PetscInfo2(snes,"Converged due to function norm %G < %G\n",fnorm,snes->abstol);CHKERRQ(ierr); 143 *reason = SNES_CONVERGED_FNORM_ABS; 144 } else if (snes->nfuncs >= snes->max_funcs) { 145 ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 146 *reason = SNES_DIVERGED_FUNCTION_COUNT; 147 } 148 149 if (it && !*reason) { 150 if (fnorm < snes->ttol) { 151 ierr = PetscInfo2(snes,"Converged due to function norm %G < %G (relative tolerance)\n",fnorm,snes->ttol);CHKERRQ(ierr); 152 *reason = SNES_CONVERGED_FNORM_RELATIVE; 153 } 154 } 155 PetscFunctionReturn(0); 156 } 157 158 /* 159 SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 160 161 Input Parameter: 162 . phi - the semismooth function 163 164 Output Parameter: 165 . merit - the merit function 166 . phinorm - ||phi|| 167 168 Notes: 169 The merit function for the mixed complementarity problem is defined as 170 merit = 0.5*phi^T*phi 171 */ 172 #undef __FUNCT__ 173 #define __FUNCT__ "SNESVIComputeMeritFunction" 174 static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm) 175 { 176 PetscErrorCode ierr; 177 178 PetscFunctionBegin; 179 ierr = VecNormBegin(phi,NORM_2,phinorm); 180 ierr = VecNormEnd(phi,NORM_2,phinorm); 181 182 *merit = 0.5*(*phinorm)*(*phinorm); 183 PetscFunctionReturn(0); 184 } 185 186 static inline PetscScalar Phi(PetscScalar a,PetscScalar b) 187 { 188 return a + b - sqrt(a*a + b*b); 189 } 190 191 static inline PetscScalar DPhi(PetscScalar a,PetscScalar b) 192 { 193 if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a/ sqrt(a*a + b*b); 194 else return .5; 195 } 196 197 /* 198 SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form. 199 200 Input Parameters: 201 . snes - the SNES context 202 . x - current iterate 203 . functx - user defined function context 204 205 Output Parameters: 206 . phi - Semismooth function 207 208 */ 209 #undef __FUNCT__ 210 #define __FUNCT__ "SNESVIComputeFunction" 211 static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx) 212 { 213 PetscErrorCode ierr; 214 SNES_VI *vi = (SNES_VI*)snes->data; 215 Vec Xl = vi->xl,Xu = vi->xu,F = snes->vec_func; 216 PetscScalar *phi_arr,*x_arr,*f_arr,*l,*u; 217 PetscInt i,nlocal; 218 219 PetscFunctionBegin; 220 ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr); 221 ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 222 ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr); 223 ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr); 224 ierr = VecGetArray(Xl,&l);CHKERRQ(ierr); 225 ierr = VecGetArray(Xu,&u);CHKERRQ(ierr); 226 ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr); 227 228 for (i=0;i < nlocal;i++) { 229 if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { /* no constraints on variable */ 230 phi_arr[i] = f_arr[i]; 231 } else if (l[i] <= PETSC_VI_NINF) { /* upper bound on variable only */ 232 phi_arr[i] = -Phi(u[i] - x_arr[i],-f_arr[i]); 233 } else if (u[i] >= PETSC_VI_INF) { /* lower bound on variable only */ 234 phi_arr[i] = Phi(x_arr[i] - l[i],f_arr[i]); 235 } else if (l[i] == u[i]) { 236 phi_arr[i] = l[i] - x_arr[i]; 237 } else { /* both bounds on variable */ 238 phi_arr[i] = Phi(x_arr[i] - l[i],-Phi(u[i] - x_arr[i],-f_arr[i])); 239 } 240 } 241 242 ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr); 243 ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr); 244 ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr); 245 ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr); 246 ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr); 247 PetscFunctionReturn(0); 248 } 249 250 /* 251 SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 252 the semismooth jacobian. 253 */ 254 #undef __FUNCT__ 255 #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors" 256 PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db) 257 { 258 PetscErrorCode ierr; 259 SNES_VI *vi = (SNES_VI*)snes->data; 260 PetscScalar *l,*u,*x,*f,*da,*db,da1,da2,db1,db2; 261 PetscInt i,nlocal; 262 263 PetscFunctionBegin; 264 265 ierr = VecGetArray(X,&x);CHKERRQ(ierr); 266 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 267 ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr); 268 ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr); 269 ierr = VecGetArray(Da,&da);CHKERRQ(ierr); 270 ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 271 ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 272 273 for (i=0;i< nlocal;i++) { 274 if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) {/* no constraints on variable */ 275 da[i] = 0; 276 db[i] = 1; 277 } else if (l[i] <= PETSC_VI_NINF) { /* upper bound on variable only */ 278 da[i] = DPhi(u[i] - x[i], -f[i]); 279 db[i] = DPhi(-f[i],u[i] - x[i]); 280 } else if (u[i] >= PETSC_VI_INF) { /* lower bound on variable only */ 281 da[i] = DPhi(x[i] - l[i], f[i]); 282 db[i] = DPhi(f[i],x[i] - l[i]); 283 } else if (l[i] == u[i]) { /* fixed variable */ 284 da[i] = 1; 285 db[i] = 0; 286 } else { /* upper and lower bounds on variable */ 287 da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i])); 288 db1 = DPhi(-Phi(u[i] - x[i], -f[i]),x[i] - l[i]); 289 da2 = DPhi(u[i] - x[i], -f[i]); 290 db2 = DPhi(-f[i],u[i] - x[i]); 291 da[i] = da1 + db1*da2; 292 db[i] = db1*db2; 293 } 294 } 295 296 ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 297 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 298 ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr); 299 ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr); 300 ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr); 301 ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 302 PetscFunctionReturn(0); 303 } 304 305 /* 306 SNESVIComputeJacobian - Computes the jacobian of the semismooth function.The Jacobian for the semismooth function is an element of the B-subdifferential of the Fischer-Burmeister function for complementarity problems. 307 308 Input Parameters: 309 . Da - Diagonal shift vector for the semismooth jacobian. 310 . Db - Row scaling vector for the semismooth jacobian. 311 312 Output Parameters: 313 . jac - semismooth jacobian 314 . jac_pre - optional preconditioning matrix 315 316 Notes: 317 The semismooth jacobian matrix is given by 318 jac = Da + Db*jacfun 319 where Db is the row scaling matrix stored as a vector, 320 Da is the diagonal perturbation matrix stored as a vector 321 and jacfun is the jacobian of the original nonlinear function. 322 */ 323 #undef __FUNCT__ 324 #define __FUNCT__ "SNESVIComputeJacobian" 325 PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db) 326 { 327 PetscErrorCode ierr; 328 329 /* Do row scaling and add diagonal perturbation */ 330 ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr); 331 ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr); 332 if (jac != jac_pre) { /* If jac and jac_pre are different */ 333 ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL); 334 ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr); 335 } 336 PetscFunctionReturn(0); 337 } 338 339 /* 340 SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi. 341 342 Input Parameters: 343 phi - semismooth function. 344 H - semismooth jacobian 345 346 Output Parameters: 347 dpsi - merit function gradient 348 349 Notes: 350 The merit function gradient is computed as follows 351 dpsi = H^T*phi 352 */ 353 #undef __FUNCT__ 354 #define __FUNCT__ "SNESVIComputeMeritFunctionGradient" 355 PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) 356 { 357 PetscErrorCode ierr; 358 359 PetscFunctionBegin; 360 ierr = MatMultTranspose(H,phi,dpsi); 361 PetscFunctionReturn(0); 362 } 363 364 /* -------------------------------------------------------------------------- */ 365 /* 366 SNESVIProjectOntoBounds - Projects X onto the feasible region so that Xl[i] <= X[i] <= Xu[i] for i = 1...n. 367 368 Input Parameters: 369 . SNES - nonlinear solver context 370 371 Output Parameters: 372 . X - Bound projected X 373 374 */ 375 376 #undef __FUNCT__ 377 #define __FUNCT__ "SNESVIProjectOntoBounds" 378 PetscErrorCode SNESVIProjectOntoBounds(SNES snes,Vec X) 379 { 380 PetscErrorCode ierr; 381 SNES_VI *vi = (SNES_VI*)snes->data; 382 const PetscScalar *xl,*xu; 383 PetscScalar *x; 384 PetscInt i,n; 385 386 PetscFunctionBegin; 387 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 388 ierr = VecGetArray(X,&x);CHKERRQ(ierr); 389 ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr); 390 ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr); 391 392 for(i = 0;i<n;i++) { 393 if (x[i] < xl[i]) x[i] = xl[i]; 394 else if (x[i] > xu[i]) x[i] = xu[i]; 395 } 396 ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 397 ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr); 398 ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr); 399 PetscFunctionReturn(0); 400 } 401 402 /* -------------------------------------------------------------------- 403 404 This file implements a semismooth truncated Newton method with a line search, 405 for solving a system of nonlinear equations in complementarity form, using the KSP, Vec, 406 and Mat interfaces for linear solvers, vectors, and matrices, 407 respectively. 408 409 The following basic routines are required for each nonlinear solver: 410 SNESCreate_XXX() - Creates a nonlinear solver context 411 SNESSetFromOptions_XXX() - Sets runtime options 412 SNESSolve_XXX() - Solves the nonlinear system 413 SNESDestroy_XXX() - Destroys the nonlinear solver context 414 The suffix "_XXX" denotes a particular implementation, in this case 415 we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving 416 systems of nonlinear equations with a line search (LS) method. 417 These routines are actually called via the common user interface 418 routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 419 SNESDestroy(), so the application code interface remains identical 420 for all nonlinear solvers. 421 422 Another key routine is: 423 SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 424 by setting data structures and options. The interface routine SNESSetUp() 425 is not usually called directly by the user, but instead is called by 426 SNESSolve() if necessary. 427 428 Additional basic routines are: 429 SNESView_XXX() - Prints details of runtime options that 430 have actually been used. 431 These are called by application codes via the interface routines 432 SNESView(). 433 434 The various types of solvers (preconditioners, Krylov subspace methods, 435 nonlinear solvers, timesteppers) are all organized similarly, so the 436 above description applies to these categories also. 437 438 -------------------------------------------------------------------- */ 439 /* 440 SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton 441 method using a line search. 442 443 Input Parameters: 444 . snes - the SNES context 445 446 Output Parameter: 447 . outits - number of iterations until termination 448 449 Application Interface Routine: SNESSolve() 450 451 Notes: 452 This implements essentially a semismooth Newton method with a 453 line search. By default a cubic backtracking line search 454 is employed, as described in the text "Numerical Methods for 455 Unconstrained Optimization and Nonlinear Equations" by Dennis 456 and Schnabel. 457 */ 458 #undef __FUNCT__ 459 #define __FUNCT__ "SNESSolveVI_SS" 460 PetscErrorCode SNESSolveVI_SS(SNES snes) 461 { 462 SNES_VI *vi = (SNES_VI*)snes->data; 463 PetscErrorCode ierr; 464 PetscInt maxits,i,lits; 465 PetscBool lssucceed; 466 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 467 PetscReal gnorm,xnorm=0,ynorm; 468 Vec Y,X,F,G,W; 469 KSPConvergedReason kspreason; 470 471 PetscFunctionBegin; 472 vi->computeuserfunction = snes->ops->computefunction; 473 snes->ops->computefunction = SNESVIComputeFunction; 474 475 snes->numFailures = 0; 476 snes->numLinearSolveFailures = 0; 477 snes->reason = SNES_CONVERGED_ITERATING; 478 479 maxits = snes->max_its; /* maximum number of iterations */ 480 X = snes->vec_sol; /* solution vector */ 481 F = snes->vec_func; /* residual vector */ 482 Y = snes->work[0]; /* work vectors */ 483 G = snes->work[1]; 484 W = snes->work[2]; 485 486 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 487 snes->iter = 0; 488 snes->norm = 0.0; 489 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 490 491 ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr); 492 ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 493 if (snes->domainerror) { 494 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 495 snes->ops->computefunction = vi->computeuserfunction; 496 PetscFunctionReturn(0); 497 } 498 /* Compute Merit function */ 499 ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 500 501 ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 502 ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 503 if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 504 505 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 506 snes->norm = vi->phinorm; 507 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 508 SNESLogConvHistory(snes,vi->phinorm,0); 509 SNESMonitor(snes,0,vi->phinorm); 510 511 /* set parameter for default relative tolerance convergence test */ 512 snes->ttol = vi->phinorm*snes->rtol; 513 /* test convergence */ 514 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 515 if (snes->reason) { 516 snes->ops->computefunction = vi->computeuserfunction; 517 PetscFunctionReturn(0); 518 } 519 520 for (i=0; i<maxits; i++) { 521 522 /* Call general purpose update function */ 523 if (snes->ops->update) { 524 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 525 } 526 527 /* Solve J Y = Phi, where J is the semismooth jacobian */ 528 /* Get the nonlinear function jacobian */ 529 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 530 /* Get the diagonal shift and row scaling vectors */ 531 ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 532 /* Compute the semismooth jacobian */ 533 ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 534 /* Compute the merit function gradient */ 535 ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 536 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 537 ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 538 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 539 540 if (kspreason < 0) { 541 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 542 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 543 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 544 break; 545 } 546 } 547 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 548 snes->linear_its += lits; 549 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 550 /* 551 if (vi->precheckstep) { 552 PetscBool changed_y = PETSC_FALSE; 553 ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 554 } 555 556 if (PetscLogPrintInfo){ 557 ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 558 } 559 */ 560 /* Compute a (scaled) negative update in the line search routine: 561 Y <- X - lambda*Y 562 and evaluate G = function(Y) (depends on the line search). 563 */ 564 ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 565 ynorm = 1; gnorm = vi->phinorm; 566 ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 567 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 568 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 569 if (snes->domainerror) { 570 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 571 snes->ops->computefunction = vi->computeuserfunction; 572 PetscFunctionReturn(0); 573 } 574 if (!lssucceed) { 575 if (++snes->numFailures >= snes->maxFailures) { 576 PetscBool ismin; 577 snes->reason = SNES_DIVERGED_LINE_SEARCH; 578 ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 579 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 580 break; 581 } 582 } 583 /* Update function and solution vectors */ 584 vi->phinorm = gnorm; 585 vi->merit = 0.5*vi->phinorm*vi->phinorm; 586 ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 587 ierr = VecCopy(W,X);CHKERRQ(ierr); 588 /* Monitor convergence */ 589 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 590 snes->iter = i+1; 591 snes->norm = vi->phinorm; 592 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 593 SNESLogConvHistory(snes,snes->norm,lits); 594 SNESMonitor(snes,snes->iter,snes->norm); 595 /* Test for convergence, xnorm = || X || */ 596 if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 597 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 598 if (snes->reason) break; 599 } 600 if (i == maxits) { 601 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 602 if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 603 } 604 snes->ops->computefunction = vi->computeuserfunction; 605 PetscFunctionReturn(0); 606 } 607 608 #undef __FUNCT__ 609 #define __FUNCT__ "SNESVICreateIndexSets_RS" 610 PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec Xl, Vec Xu,IS* ISact,IS* ISinact) 611 { 612 PetscErrorCode ierr; 613 PetscInt i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0; 614 PetscInt *idx_act,*idx_inact,i1=0,i2=0; 615 const PetscScalar *x,*xl,*xu,*f; 616 Vec F = snes->vec_func; 617 618 PetscFunctionBegin; 619 620 ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 621 ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr); 622 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 623 ierr = VecGetArrayRead(Xl,&xl);CHKERRQ(ierr); 624 ierr = VecGetArrayRead(Xu,&xu);CHKERRQ(ierr); 625 ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 626 /* Compute the sizes of the active and inactive sets */ 627 for (i=0; i < nlocal;i++) { 628 if (((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) nloc_isinact++; 629 else nloc_isact++; 630 } 631 ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 632 ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr); 633 634 /* Creating the indexing arrays */ 635 for(i=0; i < nlocal; i++) { 636 if (((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) idx_inact[i2++] = ilow+i; 637 else idx_act[i1++] = ilow+i; 638 } 639 640 /* Create the index sets */ 641 ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr); 642 ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr); 643 644 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 645 ierr = VecRestoreArrayRead(Xl,&xl);CHKERRQ(ierr); 646 ierr = VecRestoreArrayRead(Xu,&xu);CHKERRQ(ierr); 647 ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 648 ierr = PetscFree(idx_act);CHKERRQ(ierr); 649 ierr = PetscFree(idx_inact);CHKERRQ(ierr); 650 PetscFunctionReturn(0); 651 } 652 653 /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */ 654 #undef __FUNCT__ 655 #define __FUNCT__ "SNESVICreateSubVectors" 656 PetscErrorCode SNESVICreateSubVectors(SNES snes,PetscInt n,Vec* newv) 657 { 658 PetscErrorCode ierr; 659 Vec v; 660 661 PetscFunctionBegin; 662 ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr); 663 ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr); 664 ierr = VecSetFromOptions(v);CHKERRQ(ierr); 665 *newv = v; 666 667 PetscFunctionReturn(0); 668 } 669 670 /* Resets the snes PC and KSP when the active set sizes change */ 671 #undef __FUNCT__ 672 #define __FUNCT__ "SNESVIResetPCandKSP" 673 PetscErrorCode SNESVIResetPCandKSP(SNES snes,Mat Amat,Mat Pmat) 674 { 675 PetscErrorCode ierr; 676 KSP kspnew,snesksp; 677 PC pcnew; 678 const MatSolverPackage stype; 679 680 PetscFunctionBegin; 681 /* The active and inactive set sizes have changed so need to create a new snes->ksp object */ 682 ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr); 683 ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr); 684 /* Copy over snes->ksp info */ 685 kspnew->pc_side = snesksp->pc_side; 686 kspnew->rtol = snesksp->rtol; 687 kspnew->abstol = snesksp->abstol; 688 kspnew->max_it = snesksp->max_it; 689 ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr); 690 ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr); 691 ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr); 692 ierr = PCSetOperators(kspnew->pc,Amat,Pmat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 693 ierr = PCFactorGetMatSolverPackage(snesksp->pc,&stype);CHKERRQ(ierr); 694 ierr = PCFactorSetMatSolverPackage(kspnew->pc,stype);CHKERRQ(ierr); 695 ierr = KSPDestroy(snesksp);CHKERRQ(ierr); 696 snes->ksp = kspnew; 697 ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr); 698 ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr); 699 PetscFunctionReturn(0); 700 } 701 702 703 #undef __FUNCT__ 704 #define __FUNCT__ "SNESVIComputeInactiveSetFnorm" 705 PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES snes,Vec F,Vec X,PetscScalar *fnorm) 706 { 707 PetscErrorCode ierr; 708 SNES_VI *vi = (SNES_VI*)snes->data; 709 const PetscScalar *x,*xl,*xu,*f; 710 PetscInt i,n; 711 PetscReal rnorm; 712 713 PetscFunctionBegin; 714 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 715 ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr); 716 ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr); 717 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 718 ierr = VecGetArrayRead(F,&f);CHKERRQ(ierr); 719 rnorm = 0.0; 720 for (i=0; i<n; i++) { 721 if (((x[i] > xl[i] + 1.e-8 || (f[i] < 0.0)) && ((x[i] < xu[i] - 1.e-8) || f[i] > 0.0))) rnorm += f[i]*f[i]; 722 } 723 ierr = VecRestoreArrayRead(F,&f);CHKERRQ(ierr); 724 ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr); 725 ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr); 726 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 727 ierr = MPI_Allreduce(&rnorm,fnorm,1,MPIU_REAL,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr); 728 *fnorm = sqrt(*fnorm); 729 PetscFunctionReturn(0); 730 } 731 732 /* Variational Inequality solver using reduce space method. No semismooth algorithm is 733 implemented in this algorithm. It basically identifies the active variables and does 734 a linear solve on the inactive variables. */ 735 #undef __FUNCT__ 736 #define __FUNCT__ "SNESSolveVI_RS" 737 PetscErrorCode SNESSolveVI_RS(SNES snes) 738 { 739 SNES_VI *vi = (SNES_VI*)snes->data; 740 PetscErrorCode ierr; 741 PetscInt maxits,i,lits,Nis_inact; 742 PetscBool lssucceed; 743 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 744 PetscReal fnorm,gnorm,xnorm=0,ynorm; 745 Vec Y,X,F,G,W; 746 KSPConvergedReason kspreason; 747 748 PetscFunctionBegin; 749 snes->numFailures = 0; 750 snes->numLinearSolveFailures = 0; 751 snes->reason = SNES_CONVERGED_ITERATING; 752 753 maxits = snes->max_its; /* maximum number of iterations */ 754 X = snes->vec_sol; /* solution vector */ 755 F = snes->vec_func; /* residual vector */ 756 Y = snes->work[0]; /* work vectors */ 757 G = snes->work[1]; 758 W = snes->work[2]; 759 760 Nis_inact = F->map->N; 761 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 762 snes->iter = 0; 763 snes->norm = 0.0; 764 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 765 766 ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr); 767 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 768 if (snes->domainerror) { 769 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 770 PetscFunctionReturn(0); 771 } 772 ierr = SNESVIComputeInactiveSetFnorm(snes,F,X,&fnorm);CHKERRQ(ierr); 773 ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 774 ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 775 if PetscIsInfOrNanReal(fnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 776 777 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 778 snes->norm = fnorm; 779 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 780 SNESLogConvHistory(snes,fnorm,0); 781 SNESMonitor(snes,0,fnorm); 782 783 /* set parameter for default relative tolerance convergence test */ 784 snes->ttol = fnorm*snes->rtol; 785 /* test convergence */ 786 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 787 if (snes->reason) PetscFunctionReturn(0); 788 789 for (i=0; i<maxits; i++) { 790 791 IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 792 VecScatter scat_act,scat_inact; 793 PetscInt nis_act,nis_inact,Nis_inact_prev; 794 Vec Y_act,Y_inact,F_inact; 795 Mat jac_inact_inact,prejac_inact_inact; 796 797 /* Call general purpose update function */ 798 if (snes->ops->update) { 799 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 800 } 801 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 802 /* Create active and inactive index sets */ 803 ierr = SNESVICreateIndexSets_RS(snes,X,vi->xl,vi->xu,&IS_act,&IS_inact);CHKERRQ(ierr); 804 805 Nis_inact_prev = Nis_inact; 806 /* Get sizes of active and inactive sets */ 807 ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr); 808 ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr); 809 ierr = ISGetSize(IS_inact,&Nis_inact);CHKERRQ(ierr); 810 811 ierr = PetscPrintf(PETSC_COMM_WORLD,"Size of active set = %d, size of inactive set = %d\n",nis_act,nis_inact);CHKERRQ(ierr); 812 813 /* Create active and inactive set vectors */ 814 ierr = SNESVICreateSubVectors(snes,nis_inact,&F_inact);CHKERRQ(ierr); 815 ierr = SNESVICreateSubVectors(snes,nis_act,&Y_act);CHKERRQ(ierr); 816 ierr = SNESVICreateSubVectors(snes,nis_inact,&Y_inact);CHKERRQ(ierr); 817 818 /* Create inactive set submatrices */ 819 ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr); 820 821 /* Create scatter contexts */ 822 ierr = VecScatterCreate(Y,IS_act,Y_act,PETSC_NULL,&scat_act);CHKERRQ(ierr); 823 ierr = VecScatterCreate(Y,IS_inact,Y_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr); 824 825 /* Do a vec scatter to active and inactive set vectors */ 826 ierr = VecScatterBegin(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 827 ierr = VecScatterEnd(scat_inact,F,F_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 828 829 ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 830 ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 831 832 ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 833 ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 834 835 /* Active set direction = 0 */ 836 ierr = VecSet(Y_act,0);CHKERRQ(ierr); 837 if (snes->jacobian != snes->jacobian_pre) { 838 ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr); 839 } else prejac_inact_inact = jac_inact_inact; 840 841 if (Nis_inact != Nis_inact_prev) { 842 ierr = SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);CHKERRQ(ierr); 843 } 844 845 ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr); 846 ierr = SNES_KSPSolve(snes,snes->ksp,F_inact,Y_inact);CHKERRQ(ierr); 847 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 848 if (kspreason < 0) { 849 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 850 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 851 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 852 break; 853 } 854 } 855 856 ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 857 ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 858 ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 859 ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 860 861 ierr = VecDestroy(F_inact);CHKERRQ(ierr); 862 ierr = VecDestroy(Y_act);CHKERRQ(ierr); 863 ierr = VecDestroy(Y_inact);CHKERRQ(ierr); 864 ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr); 865 ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr); 866 ierr = ISDestroy(IS_act);CHKERRQ(ierr); 867 ierr = ISDestroy(IS_inact);CHKERRQ(ierr); 868 ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr); 869 if (snes->jacobian != snes->jacobian_pre) { 870 ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr); 871 } 872 873 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 874 snes->linear_its += lits; 875 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 876 /* 877 if (vi->precheckstep) { 878 PetscBool changed_y = PETSC_FALSE; 879 ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 880 } 881 882 if (PetscLogPrintInfo){ 883 ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 884 } 885 */ 886 /* Compute a (scaled) negative update in the line search routine: 887 Y <- X - lambda*Y 888 and evaluate G = function(Y) (depends on the line search). 889 */ 890 ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 891 ynorm = 1; gnorm = fnorm; 892 ierr = (*vi->LineSearch)(snes,vi->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 893 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 894 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 895 if (snes->domainerror) { 896 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 897 PetscFunctionReturn(0); 898 } 899 if (!lssucceed) { 900 if (++snes->numFailures >= snes->maxFailures) { 901 PetscBool ismin; 902 snes->reason = SNES_DIVERGED_LINE_SEARCH; 903 ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 904 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 905 break; 906 } 907 } 908 /* Update function and solution vectors */ 909 fnorm = gnorm; 910 ierr = VecCopy(G,F);CHKERRQ(ierr); 911 ierr = VecCopy(W,X);CHKERRQ(ierr); 912 /* Monitor convergence */ 913 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 914 snes->iter = i+1; 915 snes->norm = fnorm; 916 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 917 SNESLogConvHistory(snes,snes->norm,lits); 918 SNESMonitor(snes,snes->iter,snes->norm); 919 /* Test for convergence, xnorm = || X || */ 920 if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 921 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 922 if (snes->reason) break; 923 } 924 if (i == maxits) { 925 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 926 if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 927 } 928 PetscFunctionReturn(0); 929 } 930 931 /* -------------------------------------------------------------------------- */ 932 /* 933 SNESSetUp_VI - Sets up the internal data structures for the later use 934 of the SNESVI nonlinear solver. 935 936 Input Parameter: 937 . snes - the SNES context 938 . x - the solution vector 939 940 Application Interface Routine: SNESSetUp() 941 942 Notes: 943 For basic use of the SNES solvers, the user need not explicitly call 944 SNESSetUp(), since these actions will automatically occur during 945 the call to SNESSolve(). 946 */ 947 #undef __FUNCT__ 948 #define __FUNCT__ "SNESSetUp_VI" 949 PetscErrorCode SNESSetUp_VI(SNES snes) 950 { 951 PetscErrorCode ierr; 952 SNES_VI *vi = (SNES_VI*) snes->data; 953 PetscInt i_start[3],i_end[3]; 954 955 PetscFunctionBegin; 956 if (!snes->vec_sol_update) { 957 ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 958 ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 959 } 960 if (!snes->work) { 961 snes->nwork = 3; 962 ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 963 ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 964 } 965 966 /* If the lower and upper bound on variables are not set, set it to 967 -Inf and Inf */ 968 if (!vi->xl && !vi->xu) { 969 vi->usersetxbounds = PETSC_FALSE; 970 ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr); 971 ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr); 972 ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr); 973 ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr); 974 } else { 975 /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 976 ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 977 ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr); 978 ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr); 979 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])) 980 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."); 981 } 982 983 if (snes->ops->solve != SNESSolveVI_RS) { 984 ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr); 985 ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr); 986 ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr); 987 ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr); 988 ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 989 ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr); 990 991 } 992 PetscFunctionReturn(0); 993 } 994 /* -------------------------------------------------------------------------- */ 995 /* 996 SNESDestroy_VI - Destroys the private SNES_VI context that was created 997 with SNESCreate_VI(). 998 999 Input Parameter: 1000 . snes - the SNES context 1001 1002 Application Interface Routine: SNESDestroy() 1003 */ 1004 #undef __FUNCT__ 1005 #define __FUNCT__ "SNESDestroy_VI" 1006 PetscErrorCode SNESDestroy_VI(SNES snes) 1007 { 1008 SNES_VI *vi = (SNES_VI*) snes->data; 1009 PetscErrorCode ierr; 1010 1011 PetscFunctionBegin; 1012 if (snes->vec_sol_update) { 1013 ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 1014 snes->vec_sol_update = PETSC_NULL; 1015 } 1016 if (snes->nwork) { 1017 ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 1018 snes->nwork = 0; 1019 snes->work = PETSC_NULL; 1020 } 1021 1022 if (snes->ops->solve != SNESSolveVI_RS) { 1023 /* clear vectors */ 1024 ierr = VecDestroy(vi->dpsi);CHKERRQ(ierr); 1025 ierr = VecDestroy(vi->phi); CHKERRQ(ierr); 1026 ierr = VecDestroy(vi->Da); CHKERRQ(ierr); 1027 ierr = VecDestroy(vi->Db); CHKERRQ(ierr); 1028 ierr = VecDestroy(vi->z); CHKERRQ(ierr); 1029 ierr = VecDestroy(vi->t); CHKERRQ(ierr); 1030 } 1031 1032 if (!vi->usersetxbounds) { 1033 ierr = VecDestroy(vi->xl); CHKERRQ(ierr); 1034 ierr = VecDestroy(vi->xu); CHKERRQ(ierr); 1035 } 1036 1037 if (vi->lsmonitor) { 1038 ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 1039 } 1040 ierr = PetscFree(snes->data);CHKERRQ(ierr); 1041 1042 /* clear composed functions */ 1043 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 1044 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 1045 PetscFunctionReturn(0); 1046 } 1047 1048 /* -------------------------------------------------------------------------- */ 1049 #undef __FUNCT__ 1050 #define __FUNCT__ "SNESLineSearchNo_VI" 1051 1052 /* 1053 This routine does not actually do a line search but it takes a full newton 1054 step while ensuring that the new iterates remain within the constraints. 1055 1056 */ 1057 PetscErrorCode SNESLineSearchNo_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 1058 { 1059 PetscErrorCode ierr; 1060 SNES_VI *vi = (SNES_VI*)snes->data; 1061 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 1062 1063 PetscFunctionBegin; 1064 *flag = PETSC_TRUE; 1065 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1066 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 1067 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 1068 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1069 if (vi->postcheckstep) { 1070 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1071 } 1072 if (changed_y) { 1073 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 1074 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1075 } 1076 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1077 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1078 if (!snes->domainerror) { 1079 if (snes->ops->solve == SNESSolveVI_RS) { 1080 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 1081 } else { 1082 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 1083 } 1084 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1085 } 1086 if (vi->lsmonitor) { 1087 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1088 } 1089 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1090 PetscFunctionReturn(0); 1091 } 1092 1093 /* -------------------------------------------------------------------------- */ 1094 #undef __FUNCT__ 1095 #define __FUNCT__ "SNESLineSearchNoNorms_VI" 1096 1097 /* 1098 This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 1099 */ 1100 PetscErrorCode SNESLineSearchNoNorms_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 1101 { 1102 PetscErrorCode ierr; 1103 SNES_VI *vi = (SNES_VI*)snes->data; 1104 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 1105 1106 PetscFunctionBegin; 1107 *flag = PETSC_TRUE; 1108 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1109 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 1110 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1111 if (vi->postcheckstep) { 1112 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1113 } 1114 if (changed_y) { 1115 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 1116 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1117 } 1118 1119 /* don't evaluate function the last time through */ 1120 if (snes->iter < snes->max_its-1) { 1121 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1122 } 1123 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1124 PetscFunctionReturn(0); 1125 } 1126 1127 /* -------------------------------------------------------------------------- */ 1128 #undef __FUNCT__ 1129 #define __FUNCT__ "SNESLineSearchCubic_VI" 1130 /* 1131 This routine implements a cubic line search while doing a projection on the variable bounds 1132 */ 1133 PetscErrorCode SNESLineSearchCubic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 1134 { 1135 PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 1136 PetscReal minlambda,lambda,lambdatemp; 1137 #if defined(PETSC_USE_COMPLEX) 1138 PetscScalar cinitslope; 1139 #endif 1140 PetscErrorCode ierr; 1141 PetscInt count; 1142 SNES_VI *vi = (SNES_VI*)snes->data; 1143 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 1144 MPI_Comm comm; 1145 1146 PetscFunctionBegin; 1147 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 1148 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1149 *flag = PETSC_TRUE; 1150 1151 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 1152 if (*ynorm == 0.0) { 1153 if (vi->lsmonitor) { 1154 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 1155 } 1156 *gnorm = fnorm; 1157 ierr = VecCopy(x,w);CHKERRQ(ierr); 1158 ierr = VecCopy(f,g);CHKERRQ(ierr); 1159 *flag = PETSC_FALSE; 1160 goto theend1; 1161 } 1162 if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1163 if (vi->lsmonitor) { 1164 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 1165 } 1166 ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 1167 *ynorm = vi->maxstep; 1168 } 1169 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 1170 minlambda = vi->minlambda/rellength; 1171 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 1172 #if defined(PETSC_USE_COMPLEX) 1173 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 1174 initslope = PetscRealPart(cinitslope); 1175 #else 1176 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 1177 #endif 1178 if (initslope > 0.0) initslope = -initslope; 1179 if (initslope == 0.0) initslope = -1.0; 1180 1181 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1182 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1183 if (snes->nfuncs >= snes->max_funcs) { 1184 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 1185 *flag = PETSC_FALSE; 1186 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1187 goto theend1; 1188 } 1189 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1190 if (snes->ops->solve == SNESSolveVI_RS) { 1191 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 1192 } else { 1193 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1194 } 1195 if (snes->domainerror) { 1196 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1197 PetscFunctionReturn(0); 1198 } 1199 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1200 ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1201 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1202 if (vi->lsmonitor) { 1203 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1204 } 1205 goto theend1; 1206 } 1207 1208 /* Fit points with quadratic */ 1209 lambda = 1.0; 1210 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 1211 lambdaprev = lambda; 1212 gnormprev = *gnorm; 1213 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 1214 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 1215 else lambda = lambdatemp; 1216 1217 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 1218 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1219 if (snes->nfuncs >= snes->max_funcs) { 1220 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 1221 *flag = PETSC_FALSE; 1222 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1223 goto theend1; 1224 } 1225 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1226 if (snes->ops->solve == SNESSolveVI_RS) { 1227 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 1228 } else { 1229 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1230 } 1231 if (snes->domainerror) { 1232 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1233 PetscFunctionReturn(0); 1234 } 1235 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1236 if (vi->lsmonitor) { 1237 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 1238 } 1239 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1240 if (vi->lsmonitor) { 1241 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 1242 } 1243 goto theend1; 1244 } 1245 1246 /* Fit points with cubic */ 1247 count = 1; 1248 while (PETSC_TRUE) { 1249 if (lambda <= minlambda) { 1250 if (vi->lsmonitor) { 1251 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 1252 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,minlambda,lambda,initslope);CHKERRQ(ierr); 1253 } 1254 *flag = PETSC_FALSE; 1255 break; 1256 } 1257 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 1258 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 1259 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 1260 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 1261 d = b*b - 3*a*initslope; 1262 if (d < 0.0) d = 0.0; 1263 if (a == 0.0) { 1264 lambdatemp = -initslope/(2.0*b); 1265 } else { 1266 lambdatemp = (-b + sqrt(d))/(3.0*a); 1267 } 1268 lambdaprev = lambda; 1269 gnormprev = *gnorm; 1270 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 1271 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 1272 else lambda = lambdatemp; 1273 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 1274 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1275 if (snes->nfuncs >= snes->max_funcs) { 1276 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 1277 ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 1278 *flag = PETSC_FALSE; 1279 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1280 break; 1281 } 1282 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1283 if (snes->ops->solve == SNESSolveVI_RS) { 1284 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 1285 } else { 1286 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1287 } 1288 if (snes->domainerror) { 1289 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1290 PetscFunctionReturn(0); 1291 } 1292 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1293 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */ 1294 if (vi->lsmonitor) { 1295 ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 1296 } 1297 break; 1298 } else { 1299 if (vi->lsmonitor) { 1300 ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 1301 } 1302 } 1303 count++; 1304 } 1305 theend1: 1306 /* Optional user-defined check for line search step validity */ 1307 if (vi->postcheckstep && *flag) { 1308 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1309 if (changed_y) { 1310 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1311 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1312 } 1313 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 1314 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1315 if (snes->ops->solve == SNESSolveVI_RS) { 1316 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 1317 } else { 1318 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1319 } 1320 if (snes->domainerror) { 1321 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1322 PetscFunctionReturn(0); 1323 } 1324 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1325 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 1326 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 1327 } 1328 } 1329 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1330 PetscFunctionReturn(0); 1331 } 1332 1333 #undef __FUNCT__ 1334 #define __FUNCT__ "SNESLineSearchQuadratic_VI" 1335 /* 1336 This routine does a quadratic line search while keeping the iterates within the variable bounds 1337 */ 1338 PetscErrorCode SNESLineSearchQuadratic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 1339 { 1340 /* 1341 Note that for line search purposes we work with with the related 1342 minimization problem: 1343 min z(x): R^n -> R, 1344 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 1345 */ 1346 PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 1347 #if defined(PETSC_USE_COMPLEX) 1348 PetscScalar cinitslope; 1349 #endif 1350 PetscErrorCode ierr; 1351 PetscInt count; 1352 SNES_VI *vi = (SNES_VI*)snes->data; 1353 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 1354 1355 PetscFunctionBegin; 1356 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1357 *flag = PETSC_TRUE; 1358 1359 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 1360 if (*ynorm == 0.0) { 1361 if (vi->lsmonitor) { 1362 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 1363 } 1364 *gnorm = fnorm; 1365 ierr = VecCopy(x,w);CHKERRQ(ierr); 1366 ierr = VecCopy(f,g);CHKERRQ(ierr); 1367 *flag = PETSC_FALSE; 1368 goto theend2; 1369 } 1370 if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1371 ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 1372 *ynorm = vi->maxstep; 1373 } 1374 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 1375 minlambda = vi->minlambda/rellength; 1376 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 1377 #if defined(PETSC_USE_COMPLEX) 1378 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 1379 initslope = PetscRealPart(cinitslope); 1380 #else 1381 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 1382 #endif 1383 if (initslope > 0.0) initslope = -initslope; 1384 if (initslope == 0.0) initslope = -1.0; 1385 ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 1386 1387 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1388 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1389 if (snes->nfuncs >= snes->max_funcs) { 1390 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 1391 *flag = PETSC_FALSE; 1392 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1393 goto theend2; 1394 } 1395 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1396 if (snes->ops->solve == SNESSolveVI_RS) { 1397 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 1398 } else { 1399 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1400 } 1401 if (snes->domainerror) { 1402 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1403 PetscFunctionReturn(0); 1404 } 1405 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1406 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1407 if (vi->lsmonitor) { 1408 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1409 } 1410 goto theend2; 1411 } 1412 1413 /* Fit points with quadratic */ 1414 lambda = 1.0; 1415 count = 1; 1416 while (PETSC_TRUE) { 1417 if (lambda <= minlambda) { /* bad luck; use full step */ 1418 if (vi->lsmonitor) { 1419 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 1420 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 1421 } 1422 ierr = VecCopy(x,w);CHKERRQ(ierr); 1423 *flag = PETSC_FALSE; 1424 break; 1425 } 1426 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 1427 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 1428 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 1429 else lambda = lambdatemp; 1430 1431 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 1432 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1433 if (snes->nfuncs >= snes->max_funcs) { 1434 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 1435 ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 1436 *flag = PETSC_FALSE; 1437 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1438 break; 1439 } 1440 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1441 if (snes->domainerror) { 1442 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1443 PetscFunctionReturn(0); 1444 } 1445 if (snes->ops->solve == SNESSolveVI_RS) { 1446 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 1447 } else { 1448 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1449 } 1450 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1451 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1452 if (vi->lsmonitor) { 1453 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 1454 } 1455 break; 1456 } 1457 count++; 1458 } 1459 theend2: 1460 /* Optional user-defined check for line search step validity */ 1461 if (vi->postcheckstep) { 1462 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1463 if (changed_y) { 1464 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1465 ierr = SNESVIProjectOntoBounds(snes,w);CHKERRQ(ierr); 1466 } 1467 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 1468 ierr = SNESComputeFunction(snes,w,g); 1469 if (snes->domainerror) { 1470 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1471 PetscFunctionReturn(0); 1472 } 1473 if (snes->ops->solve == SNESSolveVI_RS) { 1474 ierr = SNESVIComputeInactiveSetFnorm(snes,g,w,gnorm);CHKERRQ(ierr); 1475 } else { 1476 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1477 } 1478 1479 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 1480 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 1481 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1482 } 1483 } 1484 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1485 PetscFunctionReturn(0); 1486 } 1487 1488 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/ 1489 /* -------------------------------------------------------------------------- */ 1490 EXTERN_C_BEGIN 1491 #undef __FUNCT__ 1492 #define __FUNCT__ "SNESLineSearchSet_VI" 1493 PetscErrorCode SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx) 1494 { 1495 PetscFunctionBegin; 1496 ((SNES_VI *)(snes->data))->LineSearch = func; 1497 ((SNES_VI *)(snes->data))->lsP = lsctx; 1498 PetscFunctionReturn(0); 1499 } 1500 EXTERN_C_END 1501 1502 /* -------------------------------------------------------------------------- */ 1503 EXTERN_C_BEGIN 1504 #undef __FUNCT__ 1505 #define __FUNCT__ "SNESLineSearchSetMonitor_VI" 1506 PetscErrorCode SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg) 1507 { 1508 SNES_VI *vi = (SNES_VI*)snes->data; 1509 PetscErrorCode ierr; 1510 1511 PetscFunctionBegin; 1512 if (flg && !vi->lsmonitor) { 1513 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr); 1514 } else if (!flg && vi->lsmonitor) { 1515 ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 1516 } 1517 PetscFunctionReturn(0); 1518 } 1519 EXTERN_C_END 1520 1521 /* 1522 SNESView_VI - Prints info from the SNESVI data structure. 1523 1524 Input Parameters: 1525 . SNES - the SNES context 1526 . viewer - visualization context 1527 1528 Application Interface Routine: SNESView() 1529 */ 1530 #undef __FUNCT__ 1531 #define __FUNCT__ "SNESView_VI" 1532 static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer) 1533 { 1534 SNES_VI *vi = (SNES_VI *)snes->data; 1535 const char *cstr,*tstr; 1536 PetscErrorCode ierr; 1537 PetscBool iascii; 1538 1539 PetscFunctionBegin; 1540 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1541 if (iascii) { 1542 if (vi->LineSearch == SNESLineSearchNo_VI) cstr = "SNESLineSearchNo"; 1543 else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic"; 1544 else if (vi->LineSearch == SNESLineSearchCubic_VI) cstr = "SNESLineSearchCubic"; 1545 else cstr = "unknown"; 1546 if (snes->ops->solve == SNESSolveVI_SS) tstr = "Semismooth"; 1547 else if (snes->ops->solve == SNESSolveVI_RS) tstr = "Reduced Space"; 1548 else tstr = "unknown"; 1549 ierr = PetscViewerASCIIPrintf(viewer," VI algorithm: %s\n",tstr);CHKERRQ(ierr); 1550 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1551 ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr); 1552 } else { 1553 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name); 1554 } 1555 PetscFunctionReturn(0); 1556 } 1557 1558 #undef __FUNCT__ 1559 #define __FUNCT__ "SNESVISetVariableBounds" 1560 /*@ 1561 SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 1562 1563 Input Parameters: 1564 . snes - the SNES context. 1565 . xl - lower bound. 1566 . xu - upper bound. 1567 1568 Notes: 1569 If this routine is not called then the lower and upper bounds are set to 1570 PETSC_VI_INF and PETSC_VI_NINF respectively during SNESSetUp(). 1571 1572 @*/ 1573 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 1574 { 1575 SNES_VI *vi = (SNES_VI*)snes->data; 1576 1577 PetscFunctionBegin; 1578 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1579 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 1580 PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 1581 if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1582 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); 1583 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); 1584 1585 vi->usersetxbounds = PETSC_TRUE; 1586 vi->xl = xl; 1587 vi->xu = xu; 1588 PetscFunctionReturn(0); 1589 } 1590 /* -------------------------------------------------------------------------- */ 1591 /* 1592 SNESSetFromOptions_VI - Sets various parameters for the SNESVI method. 1593 1594 Input Parameter: 1595 . snes - the SNES context 1596 1597 Application Interface Routine: SNESSetFromOptions() 1598 */ 1599 #undef __FUNCT__ 1600 #define __FUNCT__ "SNESSetFromOptions_VI" 1601 static PetscErrorCode SNESSetFromOptions_VI(SNES snes) 1602 { 1603 SNES_VI *vi = (SNES_VI *)snes->data; 1604 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1605 const char *vies[] = {"ss","rs"}; 1606 PetscErrorCode ierr; 1607 PetscInt indx; 1608 PetscBool flg,set,flg2; 1609 1610 PetscFunctionBegin; 1611 ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 1612 ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 1613 if (flg) { 1614 ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr); 1615 } 1616 ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr); 1617 ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr); 1618 ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr); 1619 ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr); 1620 ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 1621 if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1622 ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr); 1623 if (flg2) { 1624 switch (indx) { 1625 case 0: 1626 snes->ops->solve = SNESSolveVI_SS; 1627 break; 1628 case 1: 1629 snes->ops->solve = SNESSolveVI_RS; 1630 break; 1631 } 1632 } 1633 ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"basic",&indx,&flg);CHKERRQ(ierr); 1634 if (flg) { 1635 switch (indx) { 1636 case 0: 1637 ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 1638 break; 1639 case 1: 1640 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 1641 break; 1642 case 2: 1643 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 1644 break; 1645 case 3: 1646 ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 1647 break; 1648 } 1649 } 1650 ierr = PetscOptionsTail();CHKERRQ(ierr); 1651 PetscFunctionReturn(0); 1652 } 1653 /* -------------------------------------------------------------------------- */ 1654 /*MC 1655 SNESVI - Semismooth newton method based nonlinear solver that uses a line search 1656 1657 Options Database: 1658 + -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search 1659 . -snes_vi_alpha <alpha> - Sets alpha 1660 . -snes_vi_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 1661 . -snes_vi_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1662 - -snes_vi_monitor - print information about progress of line searches 1663 1664 1665 Level: beginner 1666 1667 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 1668 SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1669 SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 1670 1671 M*/ 1672 EXTERN_C_BEGIN 1673 #undef __FUNCT__ 1674 #define __FUNCT__ "SNESCreate_VI" 1675 PetscErrorCode SNESCreate_VI(SNES snes) 1676 { 1677 PetscErrorCode ierr; 1678 SNES_VI *vi; 1679 1680 PetscFunctionBegin; 1681 snes->ops->setup = SNESSetUp_VI; 1682 snes->ops->solve = SNESSolveVI_SS; 1683 snes->ops->destroy = SNESDestroy_VI; 1684 snes->ops->setfromoptions = SNESSetFromOptions_VI; 1685 snes->ops->view = SNESView_VI; 1686 snes->ops->converged = SNESDefaultConverged_VI; 1687 1688 ierr = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr); 1689 snes->data = (void*)vi; 1690 vi->alpha = 1.e-4; 1691 vi->maxstep = 1.e8; 1692 vi->minlambda = 1.e-12; 1693 vi->LineSearch = SNESLineSearchNo_VI; 1694 vi->lsP = PETSC_NULL; 1695 vi->postcheckstep = PETSC_NULL; 1696 vi->postcheck = PETSC_NULL; 1697 vi->precheckstep = PETSC_NULL; 1698 vi->precheck = PETSC_NULL; 1699 vi->const_tol = 2.2204460492503131e-16; 1700 1701 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr); 1702 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr); 1703 1704 PetscFunctionReturn(0); 1705 } 1706 EXTERN_C_END 1707