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