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