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,PetscBool *ismin) 14 { 15 PetscReal a1; 16 PetscErrorCode ierr; 17 PetscBool 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 PetscBool 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 SNESVIComputeFunction - 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__ "SNESVIComputeFunction" 188 static PetscErrorCode SNESVIComputeFunction(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 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. 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__ "SNESVIComputeJacobian" 262 PetscErrorCode SNESVIComputeJacobian(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,PetscBool* 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 SNESSolveVI_SS - 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__ "SNESSolveVI_SS" 555 PetscErrorCode SNESSolveVI_SS(SNES snes) 556 { 557 SNES_VI *vi = (SNES_VI*)snes->data; 558 PetscErrorCode ierr; 559 PetscInt maxits,i,lits; 560 PetscBool 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 PetscBool 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 PetscBool 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 #undef __FUNCT__ 691 #define __FUNCT__ "SNESVICreateIndexSets_AS" 692 PetscErrorCode SNESVICreateIndexSets_AS(SNES snes,Vec Db,PetscScalar thresh,IS* ISact,IS* ISinact) 693 { 694 PetscErrorCode ierr; 695 PetscInt i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0; 696 PetscInt *idx_act,*idx_inact,i1=0,i2=0; 697 PetscScalar *db; 698 699 PetscFunctionBegin; 700 701 ierr = VecGetLocalSize(Db,&nlocal);CHKERRQ(ierr); 702 ierr = VecGetOwnershipRange(Db,&ilow,&ihigh);CHKERRQ(ierr); 703 ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 704 /* Compute the sizes of the active and inactive sets */ 705 for (i=0; i < nlocal;i++) { 706 if (db[i] <= thresh) nloc_isact++; 707 else nloc_isact++; 708 } 709 710 ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 711 ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr); 712 713 /* Creating the indexing arrays */ 714 for(i=ilow; i < ihigh; i++) { 715 if (db[i] <= thresh) idx_act[i1++] = i; 716 else idx_inact[i2++] = i; 717 } 718 719 /* Create the index sets */ 720 ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr); 721 ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr); 722 723 ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 724 ierr = PetscFree(idx_act);CHKERRQ(ierr); 725 ierr = PetscFree(idx_inact);CHKERRQ(ierr); 726 PetscFunctionReturn(0); 727 } 728 729 /* Variational Inequality solver using active set method */ 730 #undef __FUNCT__ 731 #define __FUNCT__ "SNESSolveVI_AS" 732 PetscErrorCode SNESSolveVI_AS(SNES snes) 733 { 734 SNES_VI *vi = (SNES_VI*)snes->data; 735 PetscErrorCode ierr; 736 PetscInt maxits,i,lits; 737 PetscBool lssucceed,changedir; 738 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 739 PetscReal gnorm,xnorm=0,ynorm; 740 Vec Y,X,F,G,W; 741 KSPConvergedReason kspreason; 742 IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 743 PetscScalar thresh,J_norm1; 744 745 PetscFunctionBegin; 746 snes->numFailures = 0; 747 snes->numLinearSolveFailures = 0; 748 snes->reason = SNES_CONVERGED_ITERATING; 749 750 maxits = snes->max_its; /* maximum number of iterations */ 751 X = snes->vec_sol; /* solution vector */ 752 F = snes->vec_func; /* residual vector */ 753 Y = snes->work[0]; /* work vectors */ 754 G = snes->work[1]; 755 W = snes->work[2]; 756 757 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 758 snes->iter = 0; 759 snes->norm = 0.0; 760 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 761 762 ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 763 ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 764 if (snes->domainerror) { 765 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 766 PetscFunctionReturn(0); 767 } 768 /* Compute Merit function */ 769 ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 770 771 ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 772 ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 773 if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 774 775 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 776 snes->norm = vi->phinorm; 777 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 778 SNESLogConvHistory(snes,vi->phinorm,0); 779 SNESMonitor(snes,0,vi->phinorm); 780 781 /* set parameter for default relative tolerance convergence test */ 782 snes->ttol = vi->phinorm*snes->rtol; 783 /* test convergence */ 784 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 785 if (snes->reason) PetscFunctionReturn(0); 786 787 for (i=0; i<maxits; i++) { 788 789 /* Call general purpose update function */ 790 if (snes->ops->update) { 791 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 792 } 793 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 794 /* Compute the threshold value for creating active and inactive sets */ 795 ierr = MatNorm(snes->jacobian,NORM_1,&J_norm1);CHKERRQ(ierr); 796 thresh = PetscMin(vi->merit,1e-2)/(1+J_norm1); 797 /* Create active and inactive index sets */ 798 ierr = SNESVICreateIndexSets_AS(snes,vi->Db,thresh,&IS_act,&IS_inact);CHKERRQ(ierr); 799 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Active set semismooth algorithm not implemented yet"); 800 ierr = VecView(vi->Db,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 801 ierr = ISView(IS_act,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 802 ierr = ISView(IS_inact,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 803 804 ierr = ISDestroy(IS_act);CHKERRQ(ierr); 805 ierr = ISDestroy(IS_inact);CHKERRQ(ierr); 806 807 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 808 ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 809 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 810 ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 811 ierr = SNESVICheckDescentDirection(snes,vi->dpsi,Y,&changedir);CHKERRQ(ierr); 812 if (kspreason < 0 || changedir) { 813 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 814 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 815 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 816 break; 817 } 818 ierr = VecCopy(vi->dpsi,Y);CHKERRQ(ierr); 819 } 820 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 821 snes->linear_its += lits; 822 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 823 /* 824 if (vi->precheckstep) { 825 PetscBool changed_y = PETSC_FALSE; 826 ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 827 } 828 829 if (PetscLogPrintInfo){ 830 ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 831 } 832 */ 833 /* Compute a (scaled) negative update in the line search routine: 834 Y <- X - lambda*Y 835 and evaluate G = function(Y) (depends on the line search). 836 */ 837 ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 838 ynorm = 1; gnorm = vi->phinorm; 839 ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 840 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 841 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 842 if (snes->domainerror) { 843 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 844 PetscFunctionReturn(0); 845 } 846 if (!lssucceed) { 847 if (++snes->numFailures >= snes->maxFailures) { 848 PetscBool ismin; 849 snes->reason = SNES_DIVERGED_LINE_SEARCH; 850 ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 851 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 852 break; 853 } 854 } 855 /* Update function and solution vectors */ 856 vi->phinorm = gnorm; 857 vi->merit = 0.5*vi->phinorm*vi->phinorm; 858 ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 859 ierr = VecCopy(W,X);CHKERRQ(ierr); 860 /* Monitor convergence */ 861 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 862 snes->iter = i+1; 863 snes->norm = vi->phinorm; 864 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 865 SNESLogConvHistory(snes,snes->norm,lits); 866 SNESMonitor(snes,snes->iter,snes->norm); 867 /* Test for convergence, xnorm = || X || */ 868 if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 869 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 870 if (snes->reason) break; 871 } 872 if (i == maxits) { 873 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 874 if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 875 } 876 PetscFunctionReturn(0); 877 } 878 879 /* -------------------------------------------------------------------------- */ 880 /* 881 SNESSetUp_VI - Sets up the internal data structures for the later use 882 of the SNESVI nonlinear solver. 883 884 Input Parameter: 885 . snes - the SNES context 886 . x - the solution vector 887 888 Application Interface Routine: SNESSetUp() 889 890 Notes: 891 For basic use of the SNES solvers, the user need not explicitly call 892 SNESSetUp(), since these actions will automatically occur during 893 the call to SNESSolve(). 894 */ 895 #undef __FUNCT__ 896 #define __FUNCT__ "SNESSetUp_VI" 897 PetscErrorCode SNESSetUp_VI(SNES snes) 898 { 899 PetscErrorCode ierr; 900 SNES_VI *vi = (SNES_VI*) snes->data; 901 PetscInt i_start[3],i_end[3]; 902 903 PetscFunctionBegin; 904 if (!snes->vec_sol_update) { 905 ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 906 ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 907 } 908 if (!snes->work) { 909 snes->nwork = 3; 910 ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 911 ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 912 } 913 914 ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr); 915 ierr = VecDuplicate(snes->vec_sol, &vi->dpsi); CHKERRQ(ierr); 916 ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr); 917 ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr); 918 ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 919 ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr); 920 921 /* If the lower and upper bound on variables are not set, set it to 922 -Inf and Inf */ 923 if (!vi->xl && !vi->xu) { 924 vi->usersetxbounds = PETSC_FALSE; 925 ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr); 926 ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr); 927 ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr); 928 ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr); 929 } else { 930 /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 931 ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 932 ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr); 933 ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr); 934 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])) 935 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."); 936 } 937 938 vi->computeuserfunction = snes->ops->computefunction; 939 vi->computeuserjacobian = snes->ops->computejacobian; 940 941 snes->ops->computefunction = SNESVIComputeFunction; 942 snes->ops->computejacobian = SNESVIComputeJacobian; 943 PetscFunctionReturn(0); 944 } 945 /* -------------------------------------------------------------------------- */ 946 /* 947 SNESDestroy_VI - Destroys the private SNES_VI context that was created 948 with SNESCreate_VI(). 949 950 Input Parameter: 951 . snes - the SNES context 952 953 Application Interface Routine: SNESDestroy() 954 */ 955 #undef __FUNCT__ 956 #define __FUNCT__ "SNESDestroy_VI" 957 PetscErrorCode SNESDestroy_VI(SNES snes) 958 { 959 SNES_VI *vi = (SNES_VI*) snes->data; 960 PetscErrorCode ierr; 961 962 PetscFunctionBegin; 963 if (snes->vec_sol_update) { 964 ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 965 snes->vec_sol_update = PETSC_NULL; 966 } 967 if (snes->nwork) { 968 ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 969 snes->nwork = 0; 970 snes->work = PETSC_NULL; 971 } 972 973 /* clear vectors */ 974 ierr = VecDestroy(vi->phi); CHKERRQ(ierr); 975 ierr = VecDestroy(vi->dpsi); CHKERRQ(ierr); 976 ierr = VecDestroy(vi->Da); CHKERRQ(ierr); 977 ierr = VecDestroy(vi->Db); CHKERRQ(ierr); 978 ierr = VecDestroy(vi->z); CHKERRQ(ierr); 979 ierr = VecDestroy(vi->t); CHKERRQ(ierr); 980 if (!vi->usersetxbounds) { 981 ierr = VecDestroy(vi->xl); CHKERRQ(ierr); 982 ierr = VecDestroy(vi->xu); CHKERRQ(ierr); 983 } 984 if (vi->lsmonitor) { 985 ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 986 } 987 ierr = PetscFree(snes->data);CHKERRQ(ierr); 988 989 /* clear composed functions */ 990 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 991 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 992 PetscFunctionReturn(0); 993 } 994 /* -------------------------------------------------------------------------- */ 995 #undef __FUNCT__ 996 #define __FUNCT__ "SNESLineSearchNo_VI" 997 998 /* 999 This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c 1000 1001 */ 1002 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) 1003 { 1004 PetscErrorCode ierr; 1005 SNES_VI *vi = (SNES_VI*)snes->data; 1006 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 1007 1008 PetscFunctionBegin; 1009 *flag = PETSC_TRUE; 1010 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1011 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 1012 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 1013 if (vi->postcheckstep) { 1014 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1015 } 1016 if (changed_y) { 1017 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 1018 } 1019 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1020 if (!snes->domainerror) { 1021 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 1022 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1023 } 1024 if (vi->lsmonitor) { 1025 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1026 } 1027 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1028 PetscFunctionReturn(0); 1029 } 1030 1031 /* -------------------------------------------------------------------------- */ 1032 #undef __FUNCT__ 1033 #define __FUNCT__ "SNESLineSearchNoNorms_VI" 1034 1035 /* 1036 This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 1037 */ 1038 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) 1039 { 1040 PetscErrorCode ierr; 1041 SNES_VI *vi = (SNES_VI*)snes->data; 1042 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 1043 1044 PetscFunctionBegin; 1045 *flag = PETSC_TRUE; 1046 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1047 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 1048 if (vi->postcheckstep) { 1049 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1050 } 1051 if (changed_y) { 1052 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 1053 } 1054 1055 /* don't evaluate function the last time through */ 1056 if (snes->iter < snes->max_its-1) { 1057 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1058 } 1059 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1060 PetscFunctionReturn(0); 1061 } 1062 /* -------------------------------------------------------------------------- */ 1063 #undef __FUNCT__ 1064 #define __FUNCT__ "SNESLineSearchCubic_VI" 1065 /* 1066 This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c 1067 */ 1068 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) 1069 { 1070 /* 1071 Note that for line search purposes we work with with the related 1072 minimization problem: 1073 min z(x): R^n -> R, 1074 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 1075 This function z(x) is same as the merit function for the complementarity problem. 1076 */ 1077 1078 PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 1079 PetscReal minlambda,lambda,lambdatemp; 1080 #if defined(PETSC_USE_COMPLEX) 1081 PetscScalar cinitslope; 1082 #endif 1083 PetscErrorCode ierr; 1084 PetscInt count; 1085 SNES_VI *vi = (SNES_VI*)snes->data; 1086 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 1087 MPI_Comm comm; 1088 1089 PetscFunctionBegin; 1090 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 1091 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1092 *flag = PETSC_TRUE; 1093 1094 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 1095 if (*ynorm == 0.0) { 1096 if (vi->lsmonitor) { 1097 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 1098 } 1099 *gnorm = fnorm; 1100 ierr = VecCopy(x,w);CHKERRQ(ierr); 1101 ierr = VecCopy(f,g);CHKERRQ(ierr); 1102 *flag = PETSC_FALSE; 1103 goto theend1; 1104 } 1105 if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1106 if (vi->lsmonitor) { 1107 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 1108 } 1109 ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 1110 *ynorm = vi->maxstep; 1111 } 1112 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 1113 minlambda = vi->minlambda/rellength; 1114 /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 1115 #if defined(PETSC_USE_COMPLEX) 1116 ierr = VecDot(vi->dpsi,y,&cinitslope);CHKERRQ(ierr); 1117 initslope = PetscRealPart(cinitslope); 1118 #else 1119 ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 1120 #endif 1121 if (initslope > 0.0) initslope = -initslope; 1122 if (initslope == 0.0) initslope = -1.0; 1123 1124 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1125 if (snes->nfuncs >= snes->max_funcs) { 1126 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 1127 *flag = PETSC_FALSE; 1128 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1129 goto theend1; 1130 } 1131 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1132 if (snes->domainerror) { 1133 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1134 PetscFunctionReturn(0); 1135 } 1136 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1137 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1138 ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1139 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1140 if (vi->lsmonitor) { 1141 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1142 } 1143 goto theend1; 1144 } 1145 1146 /* Fit points with quadratic */ 1147 lambda = 1.0; 1148 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 1149 lambdaprev = lambda; 1150 gnormprev = *gnorm; 1151 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 1152 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 1153 else lambda = lambdatemp; 1154 1155 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 1156 if (snes->nfuncs >= snes->max_funcs) { 1157 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 1158 *flag = PETSC_FALSE; 1159 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1160 goto theend1; 1161 } 1162 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1163 if (snes->domainerror) { 1164 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1165 PetscFunctionReturn(0); 1166 } 1167 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1168 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1169 if (vi->lsmonitor) { 1170 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 1171 } 1172 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1173 if (vi->lsmonitor) { 1174 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 1175 } 1176 goto theend1; 1177 } 1178 1179 /* Fit points with cubic */ 1180 count = 1; 1181 while (PETSC_TRUE) { 1182 if (lambda <= minlambda) { 1183 if (vi->lsmonitor) { 1184 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 1185 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); 1186 } 1187 *flag = PETSC_FALSE; 1188 break; 1189 } 1190 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 1191 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 1192 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 1193 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 1194 d = b*b - 3*a*initslope; 1195 if (d < 0.0) d = 0.0; 1196 if (a == 0.0) { 1197 lambdatemp = -initslope/(2.0*b); 1198 } else { 1199 lambdatemp = (-b + sqrt(d))/(3.0*a); 1200 } 1201 lambdaprev = lambda; 1202 gnormprev = *gnorm; 1203 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 1204 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 1205 else lambda = lambdatemp; 1206 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 1207 if (snes->nfuncs >= snes->max_funcs) { 1208 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 1209 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); 1210 *flag = PETSC_FALSE; 1211 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1212 break; 1213 } 1214 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1215 if (snes->domainerror) { 1216 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1217 PetscFunctionReturn(0); 1218 } 1219 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1220 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1221 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */ 1222 if (vi->lsmonitor) { 1223 ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 1224 } 1225 break; 1226 } else { 1227 if (vi->lsmonitor) { 1228 ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 1229 } 1230 } 1231 count++; 1232 } 1233 theend1: 1234 /* Optional user-defined check for line search step validity */ 1235 if (vi->postcheckstep && *flag) { 1236 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1237 if (changed_y) { 1238 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1239 } 1240 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 1241 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1242 if (snes->domainerror) { 1243 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1244 PetscFunctionReturn(0); 1245 } 1246 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 1247 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1248 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 1249 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 1250 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 1251 } 1252 } 1253 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1254 PetscFunctionReturn(0); 1255 } 1256 /* -------------------------------------------------------------------------- */ 1257 #undef __FUNCT__ 1258 #define __FUNCT__ "SNESLineSearchQuadratic_VI" 1259 /* 1260 This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c 1261 */ 1262 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) 1263 { 1264 /* 1265 Note that for line search purposes we work with with the related 1266 minimization problem: 1267 min z(x): R^n -> R, 1268 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 1269 z(x) is the same as the merit function for the complementarity problem 1270 */ 1271 PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 1272 #if defined(PETSC_USE_COMPLEX) 1273 PetscScalar cinitslope; 1274 #endif 1275 PetscErrorCode ierr; 1276 PetscInt count; 1277 SNES_VI *vi = (SNES_VI*)snes->data; 1278 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 1279 1280 PetscFunctionBegin; 1281 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1282 *flag = PETSC_TRUE; 1283 1284 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 1285 if (*ynorm == 0.0) { 1286 if (vi->lsmonitor) { 1287 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 1288 } 1289 *gnorm = fnorm; 1290 ierr = VecCopy(x,w);CHKERRQ(ierr); 1291 ierr = VecCopy(f,g);CHKERRQ(ierr); 1292 *flag = PETSC_FALSE; 1293 goto theend2; 1294 } 1295 if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1296 ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 1297 *ynorm = vi->maxstep; 1298 } 1299 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 1300 minlambda = vi->minlambda/rellength; 1301 /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 1302 #if defined(PETSC_USE_COMPLEX) 1303 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 1304 initslope = PetscRealPart(cinitslope); 1305 #else 1306 ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 1307 #endif 1308 if (initslope > 0.0) initslope = -initslope; 1309 if (initslope == 0.0) initslope = -1.0; 1310 ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 1311 1312 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1313 if (snes->nfuncs >= snes->max_funcs) { 1314 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 1315 *flag = PETSC_FALSE; 1316 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1317 goto theend2; 1318 } 1319 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1320 if (snes->domainerror) { 1321 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1322 PetscFunctionReturn(0); 1323 } 1324 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1325 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1326 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1327 if (vi->lsmonitor) { 1328 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1329 } 1330 goto theend2; 1331 } 1332 1333 /* Fit points with quadratic */ 1334 lambda = 1.0; 1335 count = 1; 1336 while (PETSC_TRUE) { 1337 if (lambda <= minlambda) { /* bad luck; use full step */ 1338 if (vi->lsmonitor) { 1339 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 1340 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); 1341 } 1342 ierr = VecCopy(x,w);CHKERRQ(ierr); 1343 *flag = PETSC_FALSE; 1344 break; 1345 } 1346 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 1347 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 1348 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 1349 else lambda = lambdatemp; 1350 1351 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 1352 if (snes->nfuncs >= snes->max_funcs) { 1353 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 1354 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); 1355 *flag = PETSC_FALSE; 1356 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1357 break; 1358 } 1359 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1360 if (snes->domainerror) { 1361 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1362 PetscFunctionReturn(0); 1363 } 1364 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1365 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1366 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1367 if (vi->lsmonitor) { 1368 ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 1369 } 1370 break; 1371 } 1372 count++; 1373 } 1374 theend2: 1375 /* Optional user-defined check for line search step validity */ 1376 if (vi->postcheckstep) { 1377 ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1378 if (changed_y) { 1379 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1380 } 1381 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 1382 ierr = SNESComputeFunction(snes,w,g); 1383 if (snes->domainerror) { 1384 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1385 PetscFunctionReturn(0); 1386 } 1387 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 1388 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 1389 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 1390 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 1391 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1392 } 1393 } 1394 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1395 PetscFunctionReturn(0); 1396 } 1397 1398 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*/ 1399 /* -------------------------------------------------------------------------- */ 1400 EXTERN_C_BEGIN 1401 #undef __FUNCT__ 1402 #define __FUNCT__ "SNESLineSearchSet_VI" 1403 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx) 1404 { 1405 PetscFunctionBegin; 1406 ((SNES_VI *)(snes->data))->LineSearch = func; 1407 ((SNES_VI *)(snes->data))->lsP = lsctx; 1408 PetscFunctionReturn(0); 1409 } 1410 EXTERN_C_END 1411 1412 /* -------------------------------------------------------------------------- */ 1413 EXTERN_C_BEGIN 1414 #undef __FUNCT__ 1415 #define __FUNCT__ "SNESLineSearchSetMonitor_VI" 1416 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg) 1417 { 1418 SNES_VI *vi = (SNES_VI*)snes->data; 1419 PetscErrorCode ierr; 1420 1421 PetscFunctionBegin; 1422 if (flg && !vi->lsmonitor) { 1423 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr); 1424 } else if (!flg && vi->lsmonitor) { 1425 ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 1426 } 1427 PetscFunctionReturn(0); 1428 } 1429 EXTERN_C_END 1430 1431 /* 1432 SNESView_VI - Prints info from the SNESVI data structure. 1433 1434 Input Parameters: 1435 . SNES - the SNES context 1436 . viewer - visualization context 1437 1438 Application Interface Routine: SNESView() 1439 */ 1440 #undef __FUNCT__ 1441 #define __FUNCT__ "SNESView_VI" 1442 static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer) 1443 { 1444 SNES_VI *vi = (SNES_VI *)snes->data; 1445 const char *cstr; 1446 PetscErrorCode ierr; 1447 PetscBool iascii; 1448 1449 PetscFunctionBegin; 1450 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1451 if (iascii) { 1452 if (vi->LineSearch == SNESLineSearchNo_VI) cstr = "SNESLineSearchNo"; 1453 else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic"; 1454 else if (vi->LineSearch == SNESLineSearchCubic_VI) cstr = "SNESLineSearchCubic"; 1455 else cstr = "unknown"; 1456 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1457 ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr); 1458 } else { 1459 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name); 1460 } 1461 PetscFunctionReturn(0); 1462 } 1463 1464 /* 1465 SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 1466 1467 Input Parameters: 1468 . snes - the SNES context. 1469 . xl - lower bound. 1470 . xu - upper bound. 1471 1472 Notes: 1473 If this routine is not called then the lower and upper bounds are set to 1474 -Infinity and Infinity respectively during SNESSetUp. 1475 */ 1476 1477 #undef __FUNCT__ 1478 #define __FUNCT__ "SNESVISetVariableBounds" 1479 PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 1480 { 1481 SNES_VI *vi = (SNES_VI*)snes->data; 1482 1483 PetscFunctionBegin; 1484 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1485 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 1486 PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 1487 1488 /* Check if SNESSetFunction is called */ 1489 if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1490 1491 /* Check if the vector sizes are compatible for lower and upper bounds */ 1492 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); 1493 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); 1494 vi->usersetxbounds = PETSC_TRUE; 1495 vi->xl = xl; 1496 vi->xu = xu; 1497 1498 PetscFunctionReturn(0); 1499 } 1500 /* -------------------------------------------------------------------------- */ 1501 /* 1502 SNESSetFromOptions_VI - Sets various parameters for the SNESVI method. 1503 1504 Input Parameter: 1505 . snes - the SNES context 1506 1507 Application Interface Routine: SNESSetFromOptions() 1508 */ 1509 #undef __FUNCT__ 1510 #define __FUNCT__ "SNESSetFromOptions_VI" 1511 static PetscErrorCode SNESSetFromOptions_VI(SNES snes) 1512 { 1513 SNES_VI *vi = (SNES_VI *)snes->data; 1514 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1515 const char *vies[] = {"ss","as"}; 1516 PetscErrorCode ierr; 1517 PetscInt indx; 1518 PetscBool flg,set,flg2; 1519 1520 PetscFunctionBegin; 1521 ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 1522 ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr); 1523 ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr); 1524 ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr); 1525 ierr = PetscOptionsReal("-snes_vi_delta","descent test fraction","None",vi->delta,&vi->delta,0);CHKERRQ(ierr); 1526 ierr = PetscOptionsReal("-snes_vi_rho","descent test power","None",vi->rho,&vi->rho,0);CHKERRQ(ierr); 1527 ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr); 1528 ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 1529 if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1530 ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr); 1531 if (flg2) { 1532 switch (indx) { 1533 case 0: 1534 snes->ops->solve = SNESSolveVI_SS; 1535 break; 1536 case 1: 1537 snes->ops->solve = SNESSolveVI_AS; 1538 break; 1539 } 1540 } 1541 ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1542 if (flg) { 1543 switch (indx) { 1544 case 0: 1545 ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 1546 break; 1547 case 1: 1548 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 1549 break; 1550 case 2: 1551 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 1552 break; 1553 case 3: 1554 ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 1555 break; 1556 } 1557 } 1558 ierr = PetscOptionsTail();CHKERRQ(ierr); 1559 PetscFunctionReturn(0); 1560 } 1561 /* -------------------------------------------------------------------------- */ 1562 /*MC 1563 SNESVI - Semismooth newton method based nonlinear solver that uses a line search 1564 1565 Options Database: 1566 + -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search 1567 . -snes_vi_alpha <alpha> - Sets alpha 1568 . -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) 1569 . -snes_vi_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1570 -snes_vi_delta <delta> - Sets the fraction used in the descent test. 1571 -snes_vi_rho <rho> - Sets the power used in the descent test. 1572 For a descent direction to be accepted it has to satisfy the condition dpsi^T*d <= -delta*||d||^rho 1573 - -snes_vi_monitor - print information about progress of line searches 1574 1575 1576 Level: beginner 1577 1578 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 1579 SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1580 SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 1581 1582 M*/ 1583 EXTERN_C_BEGIN 1584 #undef __FUNCT__ 1585 #define __FUNCT__ "SNESCreate_VI" 1586 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes) 1587 { 1588 PetscErrorCode ierr; 1589 SNES_VI *vi; 1590 1591 PetscFunctionBegin; 1592 snes->ops->setup = SNESSetUp_VI; 1593 snes->ops->solve = SNESSolveVI_SS; 1594 snes->ops->destroy = SNESDestroy_VI; 1595 snes->ops->setfromoptions = SNESSetFromOptions_VI; 1596 snes->ops->view = SNESView_VI; 1597 snes->ops->converged = SNESDefaultConverged_VI; 1598 1599 ierr = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr); 1600 snes->data = (void*)vi; 1601 vi->alpha = 1.e-4; 1602 vi->maxstep = 1.e8; 1603 vi->minlambda = 1.e-12; 1604 vi->LineSearch = SNESLineSearchCubic_VI; 1605 vi->lsP = PETSC_NULL; 1606 vi->postcheckstep = PETSC_NULL; 1607 vi->postcheck = PETSC_NULL; 1608 vi->precheckstep = PETSC_NULL; 1609 vi->precheck = PETSC_NULL; 1610 vi->rho = 2.1; 1611 vi->delta = 1e-10; 1612 vi->const_tol = 2.2204460492503131e-16; 1613 vi->computessfunction = ComputeFischerFunction; 1614 1615 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr); 1616 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr); 1617 1618 PetscFunctionReturn(0); 1619 } 1620 EXTERN_C_END 1621