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