1 #define PETSCSNES_DLL 2 3 #include "../src/snes/impls/vi/viimpl.h" 4 5 /* 6 Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 7 || 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 8 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 9 for this trick. One assumes that the probability that W is in the null space of J is very, very small. 10 */ 11 #undef __FUNCT__ 12 #define __FUNCT__ "SNESVICheckLocalMin_Private" 13 PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 14 { 15 PetscReal a1; 16 PetscErrorCode ierr; 17 PetscTruth hastranspose; 18 19 PetscFunctionBegin; 20 *ismin = PETSC_FALSE; 21 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 22 if (hastranspose) { 23 /* Compute || J^T F|| */ 24 ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 25 ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 26 ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr); 27 if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 28 } else { 29 Vec work; 30 PetscScalar result; 31 PetscReal wnorm; 32 33 ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 34 ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 35 ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 36 ierr = MatMult(A,W,work);CHKERRQ(ierr); 37 ierr = VecDot(F,work,&result);CHKERRQ(ierr); 38 ierr = VecDestroy(work);CHKERRQ(ierr); 39 a1 = PetscAbsScalar(result)/(fnorm*wnorm); 40 ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr); 41 if (a1 < 1.e-4) *ismin = PETSC_TRUE; 42 } 43 PetscFunctionReturn(0); 44 } 45 46 /* 47 Checks if J^T(F - J*X) = 0 48 */ 49 #undef __FUNCT__ 50 #define __FUNCT__ "SNESVICheckResidual_Private" 51 PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 52 { 53 PetscReal a1,a2; 54 PetscErrorCode ierr; 55 PetscTruth hastranspose; 56 57 PetscFunctionBegin; 58 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 59 if (hastranspose) { 60 ierr = MatMult(A,X,W1);CHKERRQ(ierr); 61 ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 62 63 /* Compute || J^T W|| */ 64 ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 65 ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 66 ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 67 if (a1 != 0.0) { 68 ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 69 } 70 } 71 PetscFunctionReturn(0); 72 } 73 74 /* 75 SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm. 76 77 Notes: 78 The convergence criterion currently implemented is 79 merit < abstol 80 merit < rtol*merit_initial 81 */ 82 #undef __FUNCT__ 83 #define __FUNCT__ "SNESDefaultConverged_VI" 84 PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal merit,SNESConvergedReason *reason,void *dummy) 85 { 86 PetscErrorCode ierr; 87 88 PetscFunctionBegin; 89 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 90 PetscValidPointer(reason,6); 91 92 *reason = SNES_CONVERGED_ITERATING; 93 94 if (!it) { 95 /* set parameter for default relative tolerance convergence test */ 96 snes->ttol = merit*snes->rtol; 97 } 98 if (merit != merit) { 99 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 100 *reason = SNES_DIVERGED_FNORM_NAN; 101 } else if (merit < snes->abstol) { 102 ierr = PetscInfo2(snes,"Converged due to merit function %G < %G\n",merit,snes->abstol);CHKERRQ(ierr); 103 *reason = SNES_CONVERGED_FNORM_ABS; 104 } else if (snes->nfuncs >= snes->max_funcs) { 105 ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 106 *reason = SNES_DIVERGED_FUNCTION_COUNT; 107 } 108 109 if (it && !*reason) { 110 if (merit < snes->ttol) { 111 ierr = PetscInfo2(snes,"Converged due to merit function %G < %G (relative tolerance)\n",merit,snes->ttol);CHKERRQ(ierr); 112 *reason = SNES_CONVERGED_FNORM_RELATIVE; 113 } 114 } 115 PetscFunctionReturn(0); 116 } 117 118 /* 119 SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 120 121 Input Parameter: 122 . phi - the semismooth function 123 124 Output Parameter: 125 . merit - the merit function 126 . phinorm - ||phi|| 127 128 Notes: 129 The merit function for the mixed complementarity problem is defined as 130 merit = 0.5*phi^T*phi 131 */ 132 #undef __FUNCT__ 133 #define __FUNCT__ "SNESVIComputeMeritFunction" 134 static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm) 135 { 136 PetscErrorCode ierr; 137 138 PetscFunctionBegin; 139 ierr = VecNormBegin(phi,NORM_2,phinorm); 140 ierr = VecNormEnd(phi,NORM_2,phinorm); 141 142 *merit = 0.5*(*phinorm)*(*phinorm); 143 PetscFunctionReturn(0); 144 } 145 146 /* 147 ComputeFischerFunction - Computes the semismooth fischer burmeister function for a mixed complementarity equation. 148 149 Notes: 150 The Fischer-Burmeister function is defined as 151 ff(a,b) := sqrt(a*a + b*b) - a - b 152 and is used reformulate a complementarity equation as a semismooth equation. 153 */ 154 155 #undef __FUNCT__ 156 #define __FUNCT__ "ComputeFischerFunction" 157 static PetscErrorCode ComputeFischerFunction(PetscScalar a, PetscScalar b, PetscScalar* ff) 158 { 159 PetscFunctionBegin; 160 *ff = sqrt(a*a + b*b) - a - b; 161 PetscFunctionReturn(0); 162 } 163 164 /* 165 SNESVIComputeSSFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form. 166 167 Input Parameters: 168 . snes - the SNES context 169 . x - current iterate 170 . functx - user defined function context 171 172 Output Parameters: 173 . phi - Semismooth function 174 175 Notes: 176 The result of this function is done by cases: 177 + l[i] == -infinity, u[i] == infinity -- phi[i] = -f[i] 178 . l[i] == -infinity, u[i] finite -- phi[i] = ss(u[i]-x[i], -f[i]) 179 . l[i] finite, u[i] == infinity -- phi[i] = ss(x[i]-l[i], f[i]) 180 . l[i] finite < u[i] finite -- phi[i] = phi(x[i]-l[i], ss(u[i]-x[i], -f[u])) 181 - otherwise l[i] == u[i] -- phi[i] = l[i] - x[i] 182 ss is the semismoothing function used to reformulate the nonlinear equation in complementarity 183 form to semismooth form 184 185 */ 186 #undef __FUNCT__ 187 #define __FUNCT__ "SNESVIComputeSSFunction" 188 static PetscErrorCode SNESVIComputeSSFunction(SNES snes,Vec X,Vec phi,void* functx) 189 { 190 PetscErrorCode ierr; 191 SNES_VI *vi = (SNES_VI*)snes->data; 192 Vec Xl = vi->xl,Xu = vi->xu,F = snes->vec_func; 193 PetscScalar *phi_arr,*x_arr,*f_arr,*l,*u,t; 194 PetscInt i,nlocal; 195 196 PetscFunctionBegin; 197 ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr); 198 ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 199 200 ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr); 201 ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr); 202 ierr = VecGetArray(Xl,&l);CHKERRQ(ierr); 203 ierr = VecGetArray(Xu,&u);CHKERRQ(ierr); 204 ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr); 205 206 for (i=0;i < nlocal;i++) { 207 if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 208 phi_arr[i] = -f_arr[i]; 209 } 210 else if (l[i] <= PETSC_VI_NINF) { 211 t = u[i] - x_arr[i]; 212 ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]);CHKERRQ(ierr); 213 phi_arr[i] = -phi_arr[i]; 214 } 215 else if (u[i] >= PETSC_VI_INF) { 216 t = x_arr[i] - l[i]; 217 ierr = ComputeFischerFunction(t,f_arr[i],&phi_arr[i]);CHKERRQ(ierr); 218 } 219 else if (l[i] == u[i]) { 220 phi_arr[i] = l[i] - x_arr[i]; 221 } 222 else { 223 t = u[i] - x_arr[i]; 224 ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]); 225 t = x_arr[i] - l[i]; 226 ierr = ComputeFischerFunction(t,phi_arr[i],&phi_arr[i]); 227 } 228 } 229 230 ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr); 231 ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr); 232 ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr); 233 ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr); 234 ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr); 235 236 PetscFunctionReturn(0); 237 } 238 239 /* 240 SNESVIComputeSSJacobian - 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. 241 242 Input Parameters: 243 . snes - the SNES context 244 . X - the current iterate 245 . vec_func - nonlinear function evaluated at x 246 247 Output Parameters: 248 . jac - semismooth jacobian 249 . jac_pre - optional preconditioning matrix 250 . flag - flag passed on by SNESComputeJacobian. 251 . jacctx - user provided jacobian context 252 253 Notes: 254 The semismooth jacobian matrix is given by 255 jac = Da + Db*jacfun 256 where Db is the row scaling matrix stored as a vector, 257 Da is the diagonal perturbation matrix stored as a vector 258 and jacfun is the jacobian of the original nonlinear function. 259 */ 260 #undef __FUNCT__ 261 #define __FUNCT__ "SNESVIComputeSSJacobian" 262 PetscErrorCode SNESVIComputeSSJacobian(SNES snes,Vec X,Mat *jac, Mat *jac_pre, MatStructure *flg,void* jacctx) 263 { 264 PetscErrorCode ierr; 265 SNES_VI *vi = (SNES_VI*)snes->data; 266 PetscScalar *l,*u,*x,*f,*da,*db,*z,*t,t1,t2,ci,di,ei; 267 PetscInt i,nlocal; 268 Vec F = snes->vec_func; 269 270 PetscFunctionBegin; 271 272 ierr = (*vi->computeuserjacobian)(snes,X,jac,jac_pre,flg,jacctx);CHKERRQ(ierr); 273 274 ierr = VecGetArray(X,&x);CHKERRQ(ierr); 275 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 276 ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr); 277 ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr); 278 ierr = VecGetArray(vi->Da,&da);CHKERRQ(ierr); 279 ierr = VecGetArray(vi->Db,&db);CHKERRQ(ierr); 280 ierr = VecGetArray(vi->z,&z);CHKERRQ(ierr); 281 282 ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 283 /* Set the elements of the vector z: 284 z[i] = 1 if (x[i] - l[i],f[i]) = (0,0) or (u[i] - x[i],f[i]) = (0,0) 285 else z[i] = 0 286 */ 287 for(i=0;i < nlocal;i++) { 288 da[i] = db[i] = z[i] = 0; 289 if(PetscAbsScalar(f[i]) <= vi->const_tol) { 290 if ((l[i] > PETSC_VI_NINF) && (PetscAbsScalar(x[i]-l[i]) <= vi->const_tol)) { 291 da[i] = 1; 292 z[i] = 1; 293 } 294 if ((u[i] < PETSC_VI_INF) && (PetscAbsScalar(u[i]-x[i]) <= vi->const_tol)) { 295 db[i] = 1; 296 z[i] = 1; 297 } 298 } 299 } 300 ierr = VecRestoreArray(vi->z,&z);CHKERRQ(ierr); 301 ierr = MatMult(*jac,vi->z,vi->t);CHKERRQ(ierr); 302 ierr = VecGetArray(vi->t,&t);CHKERRQ(ierr); 303 /* Compute the elements of the diagonal perturbation vector Da and row scaling vector Db */ 304 for(i=0;i< nlocal;i++) { 305 /* Free variables */ 306 if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 307 da[i] = 0; db[i] = -1; 308 } 309 /* lower bounded variables */ 310 else if (u[i] >= PETSC_VI_INF) { 311 if (da[i] >= 1) { 312 t2 = PetscScalarNorm(1,t[i]); 313 da[i] = 1/t2 - 1; 314 db[i] = t[i]/t2 - 1; 315 } else { 316 t1 = x[i] - l[i]; 317 t2 = PetscScalarNorm(t1,f[i]); 318 da[i] = t1/t2 - 1; 319 db[i] = f[i]/t2 - 1; 320 } 321 } 322 /* upper bounded variables */ 323 else if (l[i] <= PETSC_VI_NINF) { 324 if (db[i] >= 1) { 325 t2 = PetscScalarNorm(1,t[i]); 326 da[i] = -1/t2 -1; 327 db[i] = -t[i]/t2 - 1; 328 } 329 else { 330 t1 = u[i] - x[i]; 331 t2 = PetscScalarNorm(t1,f[i]); 332 da[i] = t1/t2 - 1; 333 db[i] = -f[i]/t2 - 1; 334 } 335 } 336 /* Fixed variables */ 337 else if (l[i] == u[i]) { 338 da[i] = -1; 339 db[i] = 0; 340 } 341 /* Box constrained variables */ 342 else { 343 if (db[i] >= 1) { 344 t2 = PetscScalarNorm(1,t[i]); 345 ci = 1/t2 + 1; 346 di = t[i]/t2 + 1; 347 } 348 else { 349 t1 = x[i] - u[i]; 350 t2 = PetscScalarNorm(t1,f[i]); 351 ci = t1/t2 + 1; 352 di = f[i]/t2 + 1; 353 } 354 355 if (da[i] >= 1) { 356 t1 = ci + di*t[i]; 357 t2 = PetscScalarNorm(1,t1); 358 t1 = t1/t2 - 1; 359 t2 = 1/t2 - 1; 360 } 361 else { 362 ierr = ComputeFischerFunction(u[i]-x[i],-f[i],&ei);CHKERRQ(ierr); 363 t2 = PetscScalarNorm(x[i]-l[i],ei); 364 t1 = ei/t2 - 1; 365 t2 = (x[i] - l[i])/t2 - 1; 366 } 367 368 da[i] = t2 + t1*ci; 369 db[i] = t1*di; 370 } 371 } 372 ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 373 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 374 ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr); 375 ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr); 376 ierr = VecRestoreArray(vi->Da,&da);CHKERRQ(ierr); 377 ierr = VecRestoreArray(vi->Db,&db);CHKERRQ(ierr); 378 ierr = VecRestoreArray(vi->t,&t);CHKERRQ(ierr); 379 380 /* Do row scaling and add diagonal perturbation */ 381 ierr = MatDiagonalScale(*jac,vi->Db,PETSC_NULL);CHKERRQ(ierr); 382 ierr = MatDiagonalSet(*jac,vi->Da,ADD_VALUES);CHKERRQ(ierr); 383 if (*jac != *jac_pre) { /* If jac and jac_pre are different */ 384 ierr = MatDiagonalScale(*jac_pre,vi->Db,PETSC_NULL); 385 ierr = MatDiagonalSet(*jac_pre,vi->Da,ADD_VALUES);CHKERRQ(ierr); 386 } 387 388 PetscFunctionReturn(0); 389 } 390 391 /* 392 SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi. 393 394 Input Parameters: 395 . phi - semismooth function. 396 . H - semismooth jacobian 397 398 Output Parameters: 399 . dpsi - merit function gradient 400 401 Notes: 402 The merit function gradient is computed as follows 403 dpsi = H^T*phi 404 */ 405 #undef __FUNCT__ 406 #define __FUNCT__ "SNESVIComputeMeritFunctionGradient" 407 PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) 408 { 409 PetscErrorCode ierr; 410 411 PetscFunctionBegin; 412 ierr = MatMultTranspose(H,phi,dpsi); 413 414 PetscFunctionReturn(0); 415 } 416 417 /* 418 SNESVISetDescentDirection - Sets the descent direction for the semismooth algorithm 419 420 Input Parameters: 421 . snes - the SNES context. 422 . dpsi - gradient of the merit function. 423 424 Output Parameters: 425 . flg - PETSC_TRUE if the sufficient descent condition is not satisfied. 426 427 Notes: 428 The condition for the sufficient descent direction is 429 dpsi^T*Y > delta*||Y||^rho 430 where rho, delta are as defined in the SNES_VI structure. 431 If this condition is satisfied then the existing descent direction i.e. 432 the direction given by the linear solve should be used otherwise it should be set to the negative of the merit function gradient i.e -dpsi. 433 */ 434 #undef __FUNCT__ 435 #define __FUNCT__ "SNESVICheckDescentDirection" 436 PetscErrorCode SNESVICheckDescentDirection(SNES snes,Vec dpsi, Vec Y,PetscTruth* flg) 437 { 438 PetscErrorCode ierr; 439 SNES_VI *vi = (SNES_VI*)snes->data; 440 PetscScalar dpsidotY; 441 PetscReal norm_Y,rhs; 442 const PetscReal rho = vi->rho,delta=vi->delta; 443 444 PetscFunctionBegin; 445 446 *flg = PETSC_FALSE; 447 ierr = VecDot(dpsi,Y,&dpsidotY);CHKERRQ(ierr); 448 ierr = VecNormBegin(Y,NORM_2,&norm_Y);CHKERRQ(ierr); 449 ierr = VecNormEnd(Y,NORM_2,&norm_Y);CHKERRQ(ierr); 450 451 rhs = delta*PetscPowScalar(norm_Y,rho); 452 453 if (dpsidotY <= rhs) *flg = PETSC_TRUE; 454 455 PetscFunctionReturn(0); 456 } 457 458 /* 459 SNESVIAdjustInitialGuess - Readjusts the initial guess to the SNES solver supplied by the user so that the initial guess lies inside the feasible region . 460 461 Input Parameters: 462 . lb - lower bound. 463 . ub - upper bound. 464 465 Output Parameters: 466 . X - the readjusted initial guess. 467 468 Notes: 469 The readjusted initial guess X[i] = max(max(min(l[i],X[i]),min(X[i],u[i])),min(l[i],u[i])) 470 */ 471 #undef __FUNCT__ 472 #define __FUNCT__ "SNESVIAdjustInitialGuess" 473 PetscErrorCode SNESVIAdjustInitialGuess(Vec X, Vec lb, Vec ub) 474 { 475 PetscErrorCode ierr; 476 PetscInt i,nlocal; 477 PetscScalar *x,*l,*u; 478 479 PetscFunctionBegin; 480 481 ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 482 ierr = VecGetArray(X,&x);CHKERRQ(ierr); 483 ierr = VecGetArray(lb,&l);CHKERRQ(ierr); 484 ierr = VecGetArray(ub,&u);CHKERRQ(ierr); 485 486 for(i = 0; i < nlocal; i++) { 487 x[i] = PetscMax(PetscMax(PetscMin(x[i],l[i]),PetscMin(x[i],u[i])),PetscMin(l[i],u[i])); 488 } 489 490 ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 491 ierr = VecRestoreArray(lb,&l);CHKERRQ(ierr); 492 ierr = VecRestoreArray(ub,&u);CHKERRQ(ierr); 493 494 PetscFunctionReturn(0); 495 } 496 497 /* -------------------------------------------------------------------- 498 499 This file implements a semismooth truncated Newton method with a line search, 500 for solving a system of nonlinear equations in complementarity form, using the KSP, Vec, 501 and Mat interfaces for linear solvers, vectors, and matrices, 502 respectively. 503 504 The following basic routines are required for each nonlinear solver: 505 SNESCreate_XXX() - Creates a nonlinear solver context 506 SNESSetFromOptions_XXX() - Sets runtime options 507 SNESSolve_XXX() - Solves the nonlinear system 508 SNESDestroy_XXX() - Destroys the nonlinear solver context 509 The suffix "_XXX" denotes a particular implementation, in this case 510 we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving 511 systems of nonlinear equations with a line search (LS) method. 512 These routines are actually called via the common user interface 513 routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 514 SNESDestroy(), so the application code interface remains identical 515 for all nonlinear solvers. 516 517 Another key routine is: 518 SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 519 by setting data structures and options. The interface routine SNESSetUp() 520 is not usually called directly by the user, but instead is called by 521 SNESSolve() if necessary. 522 523 Additional basic routines are: 524 SNESView_XXX() - Prints details of runtime options that 525 have actually been used. 526 These are called by application codes via the interface routines 527 SNESView(). 528 529 The various types of solvers (preconditioners, Krylov subspace methods, 530 nonlinear solvers, timesteppers) are all organized similarly, so the 531 above description applies to these categories also. 532 533 -------------------------------------------------------------------- */ 534 /* 535 SNESSolve_VI - Solves the complementarity problem with a semismooth Newton 536 method using a line search. 537 538 Input Parameters: 539 . snes - the SNES context 540 541 Output Parameter: 542 . outits - number of iterations until termination 543 544 Application Interface Routine: SNESSolve() 545 546 Notes: 547 This implements essentially a semismooth Newton method with a 548 line search. By default a cubic backtracking line search 549 is employed, as described in the text "Numerical Methods for 550 Unconstrained Optimization and Nonlinear Equations" by Dennis 551 and Schnabel. 552 */ 553 #undef __FUNCT__ 554 #define __FUNCT__ "SNESSolve_VI" 555 PetscErrorCode SNESSolve_VI(SNES snes) 556 { 557 SNES_VI *vi = (SNES_VI*)snes->data; 558 PetscErrorCode ierr; 559 PetscInt maxits,i,lits; 560 PetscTruth lssucceed,changedir; 561 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 562 PetscReal gnorm,xnorm=0,ynorm; 563 Vec Y,X,F,G,W; 564 KSPConvergedReason kspreason; 565 566 PetscFunctionBegin; 567 snes->numFailures = 0; 568 snes->numLinearSolveFailures = 0; 569 snes->reason = SNES_CONVERGED_ITERATING; 570 571 maxits = snes->max_its; /* maximum number of iterations */ 572 X = snes->vec_sol; /* solution vector */ 573 F = snes->vec_func; /* residual vector */ 574 Y = snes->work[0]; /* work vectors */ 575 G = snes->work[1]; 576 W = snes->work[2]; 577 578 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 579 snes->iter = 0; 580 snes->norm = 0.0; 581 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 582 583 ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 584 ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 585 if (snes->domainerror) { 586 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 587 PetscFunctionReturn(0); 588 } 589 /* Compute Merit function */ 590 ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 591 592 ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 593 ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 594 if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 595 596 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 597 snes->norm = vi->phinorm; 598 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 599 SNESLogConvHistory(snes,vi->phinorm,0); 600 SNESMonitor(snes,0,vi->phinorm); 601 602 /* set parameter for default relative tolerance convergence test */ 603 snes->ttol = vi->phinorm*snes->rtol; 604 /* test convergence */ 605 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 606 if (snes->reason) PetscFunctionReturn(0); 607 608 for (i=0; i<maxits; i++) { 609 610 /* Call general purpose update function */ 611 if (snes->ops->update) { 612 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 613 } 614 615 /* Solve J Y = Phi, where J is the semismooth jacobian */ 616 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 617 618 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 619 ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 620 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 621 ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 622 ierr = SNESVICheckDescentDirection(snes,vi->dpsi,Y,&changedir);CHKERRQ(ierr); 623 if (kspreason < 0 || changedir) { 624 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 625 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 626 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 627 break; 628 } 629 ierr = VecCopy(vi->dpsi,Y);CHKERRQ(ierr); 630 } 631 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 632 snes->linear_its += lits; 633 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 634 /* 635 if (vi->precheckstep) { 636 PetscTruth changed_y = PETSC_FALSE; 637 ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 638 } 639 640 if (PetscLogPrintInfo){ 641 ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 642 } 643 */ 644 /* Compute a (scaled) negative update in the line search routine: 645 Y <- X - lambda*Y 646 and evaluate G = function(Y) (depends on the line search). 647 */ 648 ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 649 ynorm = 1; gnorm = vi->phinorm; 650 ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 651 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 652 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 653 if (snes->domainerror) { 654 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 655 PetscFunctionReturn(0); 656 } 657 if (!lssucceed) { 658 if (++snes->numFailures >= snes->maxFailures) { 659 PetscTruth ismin; 660 snes->reason = SNES_DIVERGED_LINE_SEARCH; 661 ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 662 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 663 break; 664 } 665 } 666 /* Update function and solution vectors */ 667 vi->phinorm = gnorm; 668 vi->merit = 0.5*vi->phinorm*vi->phinorm; 669 ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 670 ierr = VecCopy(W,X);CHKERRQ(ierr); 671 /* Monitor convergence */ 672 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 673 snes->iter = i+1; 674 snes->norm = vi->phinorm; 675 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 676 SNESLogConvHistory(snes,snes->norm,lits); 677 SNESMonitor(snes,snes->iter,snes->norm); 678 /* Test for convergence, xnorm = || X || */ 679 if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 680 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 681 if (snes->reason) break; 682 } 683 if (i == maxits) { 684 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 685 if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 686 } 687 PetscFunctionReturn(0); 688 } 689 /* -------------------------------------------------------------------------- */ 690 /* 691 SNESSetUp_VI - Sets up the internal data structures for the later use 692 of the SNESVI nonlinear solver. 693 694 Input Parameter: 695 . snes - the SNES context 696 . x - the solution vector 697 698 Application Interface Routine: SNESSetUp() 699 700 Notes: 701 For basic use of the SNES solvers, the user need not explicitly call 702 SNESSetUp(), since these actions will automatically occur during 703 the call to SNESSolve(). 704 */ 705 #undef __FUNCT__ 706 #define __FUNCT__ "SNESSetUp_VI" 707 PetscErrorCode SNESSetUp_VI(SNES snes) 708 { 709 PetscErrorCode ierr; 710 SNES_VI *vi = (SNES_VI*) snes->data; 711 PetscInt i_start[3],i_end[3]; 712 713 PetscFunctionBegin; 714 if (!snes->vec_sol_update) { 715 ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 716 ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 717 } 718 if (!snes->work) { 719 snes->nwork = 3; 720 ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 721 ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 722 } 723 724 ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr); 725 ierr = VecDuplicate(snes->vec_sol, &vi->dpsi); CHKERRQ(ierr); 726 ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr); 727 ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr); 728 ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 729 ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr); 730 731 /* If the lower and upper bound on variables are not set, set it to 732 -Inf and Inf */ 733 if (!vi->xl && !vi->xu) { 734 vi->usersetxbounds = PETSC_FALSE; 735 ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr); 736 ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr); 737 ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr); 738 ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr); 739 } else { 740 /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 741 ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 742 ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr); 743 ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr); 744 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])) 745 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."); 746 } 747 748 vi->computeuserfunction = snes->ops->computefunction; 749 vi->computeuserjacobian = snes->ops->computejacobian; 750 751 snes->ops->computefunction = SNESVIComputeSSFunction; 752 snes->ops->computejacobian = SNESVIComputeSSJacobian; 753 PetscFunctionReturn(0); 754 } 755 /* -------------------------------------------------------------------------- */ 756 /* 757 SNESDestroy_VI - Destroys the private SNES_VI context that was created 758 with SNESCreate_VI(). 759 760 Input Parameter: 761 . snes - the SNES context 762 763 Application Interface Routine: SNESDestroy() 764 */ 765 #undef __FUNCT__ 766 #define __FUNCT__ "SNESDestroy_VI" 767 PetscErrorCode SNESDestroy_VI(SNES snes) 768 { 769 SNES_VI *vi = (SNES_VI*) snes->data; 770 PetscErrorCode ierr; 771 772 PetscFunctionBegin; 773 if (snes->vec_sol_update) { 774 ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 775 snes->vec_sol_update = PETSC_NULL; 776 } 777 if (snes->nwork) { 778 ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 779 snes->nwork = 0; 780 snes->work = PETSC_NULL; 781 } 782 783 /* clear vectors */ 784 ierr = VecDestroy(vi->phi); CHKERRQ(ierr); 785 ierr = VecDestroy(vi->dpsi); CHKERRQ(ierr); 786 ierr = VecDestroy(vi->Da); CHKERRQ(ierr); 787 ierr = VecDestroy(vi->Db); CHKERRQ(ierr); 788 ierr = VecDestroy(vi->z); CHKERRQ(ierr); 789 ierr = VecDestroy(vi->t); CHKERRQ(ierr); 790 if (!vi->usersetxbounds) { 791 ierr = VecDestroy(vi->xl); CHKERRQ(ierr); 792 ierr = VecDestroy(vi->xu); CHKERRQ(ierr); 793 } 794 if (vi->monitor) { 795 ierr = PetscViewerASCIIMonitorDestroy(vi->monitor);CHKERRQ(ierr); 796 } 797 ierr = PetscFree(snes->data);CHKERRQ(ierr); 798 799 /* clear composed functions */ 800 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 801 802 PetscFunctionReturn(0); 803 } 804 /* -------------------------------------------------------------------------- */ 805 #undef __FUNCT__ 806 #define __FUNCT__ "SNESLineSearchNo_VI" 807 808 /* 809 This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c 810 811 */ 812 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,PetscTruth *flag) 813 { 814 PetscErrorCode ierr; 815 SNES_VI *vi = (SNES_VI*)snes->data; 816 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 817 818 PetscFunctionBegin; 819 *flag = PETSC_TRUE; 820 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 821 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 822 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 823 if (vi->postcheckstep) { 824 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 825 } 826 if (changed_y) { 827 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 828 } 829 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 830 if (!snes->domainerror) { 831 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 832 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 833 } 834 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 835 PetscFunctionReturn(0); 836 } 837 838 /* -------------------------------------------------------------------------- */ 839 #undef __FUNCT__ 840 #define __FUNCT__ "SNESLineSearchNoNorms_VI" 841 842 /* 843 This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 844 */ 845 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,PetscTruth *flag) 846 { 847 PetscErrorCode ierr; 848 SNES_VI *vi = (SNES_VI*)snes->data; 849 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 850 851 PetscFunctionBegin; 852 *flag = PETSC_TRUE; 853 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 854 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 855 if (vi->postcheckstep) { 856 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 857 } 858 if (changed_y) { 859 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 860 } 861 862 /* don't evaluate function the last time through */ 863 if (snes->iter < snes->max_its-1) { 864 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 865 } 866 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 867 PetscFunctionReturn(0); 868 } 869 /* -------------------------------------------------------------------------- */ 870 #undef __FUNCT__ 871 #define __FUNCT__ "SNESLineSearchCubic_VI" 872 /* 873 This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c 874 */ 875 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,PetscTruth *flag) 876 { 877 /* 878 Note that for line search purposes we work with with the related 879 minimization problem: 880 min z(x): R^n -> R, 881 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 882 This function z(x) is same as the merit function for the complementarity problem. 883 */ 884 885 PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 886 PetscReal minlambda,lambda,lambdatemp; 887 #if defined(PETSC_USE_COMPLEX) 888 PetscScalar cinitslope; 889 #endif 890 PetscErrorCode ierr; 891 PetscInt count; 892 SNES_VI *vi = (SNES_VI*)snes->data; 893 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 894 MPI_Comm comm; 895 896 PetscFunctionBegin; 897 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 898 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 899 *flag = PETSC_TRUE; 900 901 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 902 if (*ynorm == 0.0) { 903 if (vi->monitor) { 904 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 905 } 906 *gnorm = fnorm; 907 ierr = VecCopy(x,w);CHKERRQ(ierr); 908 ierr = VecCopy(f,g);CHKERRQ(ierr); 909 *flag = PETSC_FALSE; 910 goto theend1; 911 } 912 if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 913 if (vi->monitor) { 914 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 915 } 916 ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 917 *ynorm = vi->maxstep; 918 } 919 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 920 minlambda = vi->minlambda/rellength; 921 /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 922 #if defined(PETSC_USE_COMPLEX) 923 ierr = VecDot(vi->dpsi,y,&cinitslope);CHKERRQ(ierr); 924 initslope = PetscRealPart(cinitslope); 925 #else 926 ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 927 #endif 928 if (initslope > 0.0) initslope = -initslope; 929 if (initslope == 0.0) initslope = -1.0; 930 931 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 932 if (snes->nfuncs >= snes->max_funcs) { 933 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 934 *flag = PETSC_FALSE; 935 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 936 goto theend1; 937 } 938 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 939 if (snes->domainerror) { 940 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 941 PetscFunctionReturn(0); 942 } 943 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 944 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 945 ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 946 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 947 if (vi->monitor) { 948 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 949 } 950 goto theend1; 951 } 952 953 /* Fit points with quadratic */ 954 lambda = 1.0; 955 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 956 lambdaprev = lambda; 957 gnormprev = *gnorm; 958 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 959 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 960 else lambda = lambdatemp; 961 962 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 963 if (snes->nfuncs >= snes->max_funcs) { 964 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 965 *flag = PETSC_FALSE; 966 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 967 goto theend1; 968 } 969 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 970 if (snes->domainerror) { 971 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 972 PetscFunctionReturn(0); 973 } 974 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 975 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 976 if (vi->monitor) { 977 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 978 } 979 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 980 if (vi->monitor) { 981 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 982 } 983 goto theend1; 984 } 985 986 /* Fit points with cubic */ 987 count = 1; 988 while (PETSC_TRUE) { 989 if (lambda <= minlambda) { 990 if (vi->monitor) { 991 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 992 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," 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); 993 } 994 *flag = PETSC_FALSE; 995 break; 996 } 997 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 998 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 999 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 1000 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 1001 d = b*b - 3*a*initslope; 1002 if (d < 0.0) d = 0.0; 1003 if (a == 0.0) { 1004 lambdatemp = -initslope/(2.0*b); 1005 } else { 1006 lambdatemp = (-b + sqrt(d))/(3.0*a); 1007 } 1008 lambdaprev = lambda; 1009 gnormprev = *gnorm; 1010 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 1011 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 1012 else lambda = lambdatemp; 1013 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 1014 if (snes->nfuncs >= snes->max_funcs) { 1015 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 1016 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); 1017 *flag = PETSC_FALSE; 1018 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1019 break; 1020 } 1021 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1022 if (snes->domainerror) { 1023 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1024 PetscFunctionReturn(0); 1025 } 1026 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1027 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1028 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */ 1029 if (vi->monitor) { 1030 ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 1031 } 1032 break; 1033 } else { 1034 if (vi->monitor) { 1035 ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 1036 } 1037 } 1038 count++; 1039 } 1040 theend1: 1041 /* Optional user-defined check for line search step validity */ 1042 if (vi->postcheckstep && *flag) { 1043 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1044 if (changed_y) { 1045 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1046 } 1047 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 1048 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1049 if (snes->domainerror) { 1050 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1051 PetscFunctionReturn(0); 1052 } 1053 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 1054 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1055 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 1056 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 1057 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 1058 } 1059 } 1060 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1061 PetscFunctionReturn(0); 1062 } 1063 /* -------------------------------------------------------------------------- */ 1064 #undef __FUNCT__ 1065 #define __FUNCT__ "SNESLineSearchQuadratic_VI" 1066 /* 1067 This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c 1068 */ 1069 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,PetscTruth *flag) 1070 { 1071 /* 1072 Note that for line search purposes we work with with the related 1073 minimization problem: 1074 min z(x): R^n -> R, 1075 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 1076 z(x) is the same as the merit function for the complementarity problem 1077 */ 1078 PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 1079 #if defined(PETSC_USE_COMPLEX) 1080 PetscScalar cinitslope; 1081 #endif 1082 PetscErrorCode ierr; 1083 PetscInt count; 1084 SNES_VI *vi = (SNES_VI*)snes->data; 1085 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 1086 1087 PetscFunctionBegin; 1088 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1089 *flag = PETSC_TRUE; 1090 1091 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 1092 if (*ynorm == 0.0) { 1093 if (vi->monitor) { 1094 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 1095 } 1096 *gnorm = fnorm; 1097 ierr = VecCopy(x,w);CHKERRQ(ierr); 1098 ierr = VecCopy(f,g);CHKERRQ(ierr); 1099 *flag = PETSC_FALSE; 1100 goto theend2; 1101 } 1102 if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1103 ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 1104 *ynorm = vi->maxstep; 1105 } 1106 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 1107 minlambda = vi->minlambda/rellength; 1108 /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 1109 #if defined(PETSC_USE_COMPLEX) 1110 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 1111 initslope = PetscRealPart(cinitslope); 1112 #else 1113 ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 1114 #endif 1115 if (initslope > 0.0) initslope = -initslope; 1116 if (initslope == 0.0) initslope = -1.0; 1117 ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 1118 1119 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1120 if (snes->nfuncs >= snes->max_funcs) { 1121 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 1122 *flag = PETSC_FALSE; 1123 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1124 goto theend2; 1125 } 1126 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1127 if (snes->domainerror) { 1128 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1129 PetscFunctionReturn(0); 1130 } 1131 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1132 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1133 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1134 if (vi->monitor) { 1135 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line Search: Using full step\n");CHKERRQ(ierr); 1136 } 1137 goto theend2; 1138 } 1139 1140 /* Fit points with quadratic */ 1141 lambda = 1.0; 1142 count = 1; 1143 while (PETSC_TRUE) { 1144 if (lambda <= minlambda) { /* bad luck; use full step */ 1145 if (vi->monitor) { 1146 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 1147 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor,"Line search: fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 1148 } 1149 ierr = VecCopy(x,w);CHKERRQ(ierr); 1150 *flag = PETSC_FALSE; 1151 break; 1152 } 1153 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 1154 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 1155 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 1156 else lambda = lambdatemp; 1157 1158 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 1159 if (snes->nfuncs >= snes->max_funcs) { 1160 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 1161 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); 1162 *flag = PETSC_FALSE; 1163 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1164 break; 1165 } 1166 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1167 if (snes->domainerror) { 1168 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1169 PetscFunctionReturn(0); 1170 } 1171 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1172 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1173 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1174 if (vi->monitor) { 1175 ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 1176 } 1177 break; 1178 } 1179 count++; 1180 } 1181 theend2: 1182 /* Optional user-defined check for line search step validity */ 1183 if (vi->postcheckstep) { 1184 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1185 if (changed_y) { 1186 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1187 } 1188 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 1189 ierr = SNESComputeFunction(snes,w,g); 1190 if (snes->domainerror) { 1191 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1192 PetscFunctionReturn(0); 1193 } 1194 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 1195 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 1196 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 1197 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 1198 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1199 } 1200 } 1201 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1202 PetscFunctionReturn(0); 1203 } 1204 1205 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 1206 /* -------------------------------------------------------------------------- */ 1207 EXTERN_C_BEGIN 1208 #undef __FUNCT__ 1209 #define __FUNCT__ "SNESLineSearchSet_VI" 1210 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx) 1211 { 1212 PetscFunctionBegin; 1213 ((SNES_VI *)(snes->data))->LineSearch = func; 1214 ((SNES_VI *)(snes->data))->lsP = lsctx; 1215 PetscFunctionReturn(0); 1216 } 1217 EXTERN_C_END 1218 1219 /* -------------------------------------------------------------------------- */ 1220 EXTERN_C_BEGIN 1221 #undef __FUNCT__ 1222 #define __FUNCT__ "SNESLineSearchSetMonitor_VI" 1223 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscTruth flg) 1224 { 1225 SNES_VI *vi = (SNES_VI*)snes->data; 1226 PetscErrorCode ierr; 1227 1228 PetscFunctionBegin; 1229 if (flg && !vi->monitor) { 1230 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->monitor);CHKERRQ(ierr); 1231 } else if (!flg && vi->monitor) { 1232 ierr = PetscViewerASCIIMonitorDestroy(vi->monitor);CHKERRQ(ierr); 1233 } 1234 PetscFunctionReturn(0); 1235 } 1236 EXTERN_C_END 1237 1238 /* 1239 SNESView_VI - Prints info from the SNESVI data structure. 1240 1241 Input Parameters: 1242 . SNES - the SNES context 1243 . viewer - visualization context 1244 1245 Application Interface Routine: SNESView() 1246 */ 1247 #undef __FUNCT__ 1248 #define __FUNCT__ "SNESView_VI" 1249 static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer) 1250 { 1251 SNES_VI *vi = (SNES_VI *)snes->data; 1252 const char *cstr; 1253 PetscErrorCode ierr; 1254 PetscTruth iascii; 1255 1256 PetscFunctionBegin; 1257 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1258 if (iascii) { 1259 if (vi->LineSearch == SNESLineSearchNo_VI) cstr = "SNESLineSearchNo"; 1260 else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic"; 1261 else if (vi->LineSearch == SNESLineSearchCubic_VI) cstr = "SNESLineSearchCubic"; 1262 else cstr = "unknown"; 1263 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1264 ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr); 1265 } else { 1266 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name); 1267 } 1268 PetscFunctionReturn(0); 1269 } 1270 1271 /* 1272 SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 1273 1274 Input Parameters: 1275 . snes - the SNES context. 1276 . xl - lower bound. 1277 . xu - upper bound. 1278 1279 Notes: 1280 If this routine is not called then the lower and upper bounds are set to 1281 -Infinity and Infinity respectively during SNESSetUp. 1282 */ 1283 1284 #undef __FUNCT__ 1285 #define __FUNCT__ "SNESVISetVariableBounds" 1286 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 1287 { 1288 SNES_VI *vi = (SNES_VI*)snes->data; 1289 1290 PetscFunctionBegin; 1291 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1292 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 1293 PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 1294 1295 /* Check if SNESSetFunction is called */ 1296 if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1297 1298 /* Check if the vector sizes are compatible for lower and upper bounds */ 1299 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); 1300 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); 1301 vi->usersetxbounds = PETSC_TRUE; 1302 vi->xl = xl; 1303 vi->xu = xu; 1304 1305 PetscFunctionReturn(0); 1306 } 1307 /* -------------------------------------------------------------------------- */ 1308 /* 1309 SNESSetFromOptions_VI - Sets various parameters for the SNESVI method. 1310 1311 Input Parameter: 1312 . snes - the SNES context 1313 1314 Application Interface Routine: SNESSetFromOptions() 1315 */ 1316 #undef __FUNCT__ 1317 #define __FUNCT__ "SNESSetFromOptions_VI" 1318 static PetscErrorCode SNESSetFromOptions_VI(SNES snes) 1319 { 1320 SNES_VI *vi = (SNES_VI *)snes->data; 1321 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1322 PetscErrorCode ierr; 1323 PetscInt indx; 1324 PetscTruth flg,set; 1325 1326 PetscFunctionBegin; 1327 ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 1328 ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr); 1329 ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr); 1330 ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr); 1331 ierr = PetscOptionsReal("-snes_vi_delta","descent test fraction","None",vi->delta,&vi->delta,0);CHKERRQ(ierr); 1332 ierr = PetscOptionsReal("-snes_vi_rho","descent test power","None",vi->rho,&vi->rho,0);CHKERRQ(ierr); 1333 ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr); 1334 ierr = PetscOptionsTruth("-snes_vi_monitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 1335 if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1336 1337 ierr = PetscOptionsEList("-snes_vi","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1338 if (flg) { 1339 switch (indx) { 1340 case 0: 1341 ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 1342 break; 1343 case 1: 1344 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 1345 break; 1346 case 2: 1347 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 1348 break; 1349 case 3: 1350 ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 1351 break; 1352 } 1353 } 1354 ierr = PetscOptionsTail();CHKERRQ(ierr); 1355 PetscFunctionReturn(0); 1356 } 1357 /* -------------------------------------------------------------------------- */ 1358 /*MC 1359 SNESVI - Semismooth newton method based nonlinear solver that uses a line search 1360 1361 Options Database: 1362 + -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search 1363 . -snes_vi_alpha <alpha> - Sets alpha 1364 . -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) 1365 . -snes_vi_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1366 -snes_vi_delta <delta> - Sets the fraction used in the descent test. 1367 -snes_vi_rho <rho> - Sets the power used in the descent test. 1368 For a descent direction to be accepted it has to satisfy the condition dpsi^T*d <= -delta*||d||^rho 1369 - -snes_vi_monitor - print information about progress of line searches 1370 1371 1372 Level: beginner 1373 1374 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 1375 SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1376 SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 1377 1378 M*/ 1379 EXTERN_C_BEGIN 1380 #undef __FUNCT__ 1381 #define __FUNCT__ "SNESCreate_VI" 1382 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes) 1383 { 1384 PetscErrorCode ierr; 1385 SNES_VI *vi; 1386 1387 PetscFunctionBegin; 1388 snes->ops->setup = SNESSetUp_VI; 1389 snes->ops->solve = SNESSolve_VI; 1390 snes->ops->destroy = SNESDestroy_VI; 1391 snes->ops->setfromoptions = SNESSetFromOptions_VI; 1392 snes->ops->view = SNESView_VI; 1393 snes->ops->converged = SNESDefaultConverged_VI; 1394 1395 ierr = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr); 1396 snes->data = (void*)vi; 1397 vi->alpha = 1.e-4; 1398 vi->maxstep = 1.e8; 1399 vi->minlambda = 1.e-12; 1400 vi->LineSearch = SNESLineSearchCubic_VI; 1401 vi->lsP = PETSC_NULL; 1402 vi->postcheckstep = PETSC_NULL; 1403 vi->postcheck = PETSC_NULL; 1404 vi->precheckstep = PETSC_NULL; 1405 vi->precheck = PETSC_NULL; 1406 vi->rho = 2.1; 1407 vi->delta = 1e-10; 1408 vi->const_tol = 2.2204460492503131e-16; 1409 vi->computessfunction = ComputeFischerFunction; 1410 1411 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr); 1412 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr); 1413 1414 PetscFunctionReturn(0); 1415 } 1416 EXTERN_C_END 1417