1 2 #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/ 3 4 /*@ 5 SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 6 7 Input Parameter: 8 . phi - the `Vec` holding the evaluation of the semismooth function 9 10 Output Parameters: 11 + merit - the merit function 1/2 ||phi||^2 12 - phinorm - the two-norm of the vector, ||phi|| 13 14 Level: developer 15 16 .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeFunction()` 17 @*/ 18 PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal *merit, PetscReal *phinorm) { 19 PetscFunctionBegin; 20 PetscCall(VecNormBegin(phi, NORM_2, phinorm)); 21 PetscCall(VecNormEnd(phi, NORM_2, phinorm)); 22 *merit = 0.5 * (*phinorm) * (*phinorm); 23 PetscFunctionReturn(0); 24 } 25 26 static inline PetscScalar Phi(PetscScalar a, PetscScalar b) { 27 return a + b - PetscSqrtScalar(a * a + b * b); 28 } 29 30 static inline PetscScalar DPhi(PetscScalar a, PetscScalar b) { 31 if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a / PetscSqrtScalar(a * a + b * b); 32 else return .5; 33 } 34 35 /*@ 36 SNESVIComputeFunction - Provides the function that reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear 37 equations in semismooth form. 38 39 Input Parameters: 40 + snes - the SNES context 41 . X - current iterate 42 - functx - user defined function context 43 44 Output Parameter: 45 . phi - the evaluation of Semismooth function at X 46 47 Level: developer 48 49 .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeMeritFunction()` 50 @*/ 51 PetscErrorCode SNESVIComputeFunction(SNES snes, Vec X, Vec phi, void *functx) { 52 SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 53 Vec Xl = snes->xl, Xu = snes->xu, F = snes->vec_func; 54 PetscScalar *phi_arr, *f_arr, *l, *u; 55 const PetscScalar *x_arr; 56 PetscInt i, nlocal; 57 58 PetscFunctionBegin; 59 PetscCall((*vi->computeuserfunction)(snes, X, F, functx)); 60 PetscCall(VecGetLocalSize(X, &nlocal)); 61 PetscCall(VecGetArrayRead(X, &x_arr)); 62 PetscCall(VecGetArray(F, &f_arr)); 63 PetscCall(VecGetArray(Xl, &l)); 64 PetscCall(VecGetArray(Xu, &u)); 65 PetscCall(VecGetArray(phi, &phi_arr)); 66 67 for (i = 0; i < nlocal; i++) { 68 if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */ 69 phi_arr[i] = f_arr[i]; 70 } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */ 71 phi_arr[i] = -Phi(u[i] - x_arr[i], -f_arr[i]); 72 } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */ 73 phi_arr[i] = Phi(x_arr[i] - l[i], f_arr[i]); 74 } else if (l[i] == u[i]) { 75 phi_arr[i] = l[i] - x_arr[i]; 76 } else { /* both bounds on variable */ 77 phi_arr[i] = Phi(x_arr[i] - l[i], -Phi(u[i] - x_arr[i], -f_arr[i])); 78 } 79 } 80 81 PetscCall(VecRestoreArrayRead(X, &x_arr)); 82 PetscCall(VecRestoreArray(F, &f_arr)); 83 PetscCall(VecRestoreArray(Xl, &l)); 84 PetscCall(VecRestoreArray(Xu, &u)); 85 PetscCall(VecRestoreArray(phi, &phi_arr)); 86 PetscFunctionReturn(0); 87 } 88 89 /* 90 SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 91 the semismooth jacobian. 92 */ 93 PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes, Vec X, Vec F, Mat jac, Vec Da, Vec Db) { 94 PetscScalar *l, *u, *x, *f, *da, *db, da1, da2, db1, db2; 95 PetscInt i, nlocal; 96 97 PetscFunctionBegin; 98 PetscCall(VecGetArray(X, &x)); 99 PetscCall(VecGetArray(F, &f)); 100 PetscCall(VecGetArray(snes->xl, &l)); 101 PetscCall(VecGetArray(snes->xu, &u)); 102 PetscCall(VecGetArray(Da, &da)); 103 PetscCall(VecGetArray(Db, &db)); 104 PetscCall(VecGetLocalSize(X, &nlocal)); 105 106 for (i = 0; i < nlocal; i++) { 107 if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */ 108 da[i] = 0; 109 db[i] = 1; 110 } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */ 111 da[i] = DPhi(u[i] - x[i], -f[i]); 112 db[i] = DPhi(-f[i], u[i] - x[i]); 113 } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */ 114 da[i] = DPhi(x[i] - l[i], f[i]); 115 db[i] = DPhi(f[i], x[i] - l[i]); 116 } else if (l[i] == u[i]) { /* fixed variable */ 117 da[i] = 1; 118 db[i] = 0; 119 } else { /* upper and lower bounds on variable */ 120 da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i])); 121 db1 = DPhi(-Phi(u[i] - x[i], -f[i]), x[i] - l[i]); 122 da2 = DPhi(u[i] - x[i], -f[i]); 123 db2 = DPhi(-f[i], u[i] - x[i]); 124 da[i] = da1 + db1 * da2; 125 db[i] = db1 * db2; 126 } 127 } 128 129 PetscCall(VecRestoreArray(X, &x)); 130 PetscCall(VecRestoreArray(F, &f)); 131 PetscCall(VecRestoreArray(snes->xl, &l)); 132 PetscCall(VecRestoreArray(snes->xu, &u)); 133 PetscCall(VecRestoreArray(Da, &da)); 134 PetscCall(VecRestoreArray(Db, &db)); 135 PetscFunctionReturn(0); 136 } 137 138 /* 139 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. 140 141 Input Parameters: 142 . Da - Diagonal shift vector for the semismooth jacobian. 143 . Db - Row scaling vector for the semismooth jacobian. 144 145 Output Parameters: 146 . jac - semismooth jacobian 147 . jac_pre - optional preconditioning matrix 148 149 Note: 150 The semismooth jacobian matrix is given by 151 jac = Da + Db*jacfun 152 where Db is the row scaling matrix stored as a vector, 153 Da is the diagonal perturbation matrix stored as a vector 154 and jacfun is the jacobian of the original nonlinear function. 155 */ 156 PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre, Vec Da, Vec Db) { 157 /* Do row scaling and add diagonal perturbation */ 158 PetscFunctionBegin; 159 PetscCall(MatDiagonalScale(jac, Db, NULL)); 160 PetscCall(MatDiagonalSet(jac, Da, ADD_VALUES)); 161 if (jac != jac_pre) { /* If jac and jac_pre are different */ 162 PetscCall(MatDiagonalScale(jac_pre, Db, NULL)); 163 PetscCall(MatDiagonalSet(jac_pre, Da, ADD_VALUES)); 164 } 165 PetscFunctionReturn(0); 166 } 167 168 /* 169 SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi. 170 171 Input Parameters: 172 phi - semismooth function. 173 H - semismooth jacobian 174 175 Output Parameter: 176 dpsi - merit function gradient 177 178 Note: 179 The merit function gradient is computed as follows 180 dpsi = H^T*phi 181 */ 182 PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) { 183 PetscFunctionBegin; 184 PetscCall(MatMultTranspose(H, phi, dpsi)); 185 PetscFunctionReturn(0); 186 } 187 188 /* 189 SNESSolve_VINEWTONSSLS - Solves the complementarity problem with a semismooth Newton 190 method using a line search. 191 192 Input Parameter: 193 . snes - the SNES context 194 195 Application Interface Routine: SNESSolve() 196 197 Note: 198 This implements essentially a semismooth Newton method with a 199 line search. The default line search does not do any line search 200 but rather takes a full Newton step. 201 202 Developer Note: the code in this file should be slightly modified so that this routine need not exist and the SNESSolve_NEWTONLS() routine is called directly with the appropriate wrapped function and Jacobian evaluations 203 204 */ 205 PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes) { 206 SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 207 PetscInt maxits, i, lits; 208 SNESLineSearchReason lssucceed; 209 PetscReal gnorm, xnorm = 0, ynorm; 210 Vec Y, X, F; 211 KSPConvergedReason kspreason; 212 DM dm; 213 DMSNES sdm; 214 215 PetscFunctionBegin; 216 PetscCall(SNESGetDM(snes, &dm)); 217 PetscCall(DMGetDMSNES(dm, &sdm)); 218 219 vi->computeuserfunction = sdm->ops->computefunction; 220 sdm->ops->computefunction = SNESVIComputeFunction; 221 222 snes->numFailures = 0; 223 snes->numLinearSolveFailures = 0; 224 snes->reason = SNES_CONVERGED_ITERATING; 225 226 maxits = snes->max_its; /* maximum number of iterations */ 227 X = snes->vec_sol; /* solution vector */ 228 F = snes->vec_func; /* residual vector */ 229 Y = snes->work[0]; /* work vectors */ 230 231 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 232 snes->iter = 0; 233 snes->norm = 0.0; 234 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 235 236 PetscCall(SNESVIProjectOntoBounds(snes, X)); 237 PetscCall(SNESComputeFunction(snes, X, vi->phi)); 238 if (snes->domainerror) { 239 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 240 sdm->ops->computefunction = vi->computeuserfunction; 241 PetscFunctionReturn(0); 242 } 243 /* Compute Merit function */ 244 PetscCall(SNESVIComputeMeritFunction(vi->phi, &vi->merit, &vi->phinorm)); 245 246 PetscCall(VecNormBegin(X, NORM_2, &xnorm)); /* xnorm <- ||x|| */ 247 PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 248 SNESCheckFunctionNorm(snes, vi->merit); 249 250 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 251 snes->norm = vi->phinorm; 252 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 253 PetscCall(SNESLogConvergenceHistory(snes, vi->phinorm, 0)); 254 PetscCall(SNESMonitor(snes, 0, vi->phinorm)); 255 256 /* test convergence */ 257 PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, vi->phinorm, &snes->reason, snes->cnvP); 258 if (snes->reason) { 259 sdm->ops->computefunction = vi->computeuserfunction; 260 PetscFunctionReturn(0); 261 } 262 263 for (i = 0; i < maxits; i++) { 264 /* Call general purpose update function */ 265 PetscTryTypeMethod(snes, update, snes->iter); 266 267 /* Solve J Y = Phi, where J is the semismooth jacobian */ 268 269 /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/ 270 sdm->ops->computefunction = vi->computeuserfunction; 271 PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 272 SNESCheckJacobianDomainerror(snes); 273 sdm->ops->computefunction = SNESVIComputeFunction; 274 275 /* Get the diagonal shift and row scaling vectors */ 276 PetscCall(SNESVIComputeBsubdifferentialVectors(snes, X, F, snes->jacobian, vi->Da, vi->Db)); 277 /* Compute the semismooth jacobian */ 278 PetscCall(SNESVIComputeJacobian(snes->jacobian, snes->jacobian_pre, vi->Da, vi->Db)); 279 /* Compute the merit function gradient */ 280 PetscCall(SNESVIComputeMeritFunctionGradient(snes->jacobian, vi->phi, vi->dpsi)); 281 PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 282 PetscCall(KSPSolve(snes->ksp, vi->phi, Y)); 283 PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason)); 284 285 if (kspreason < 0) { 286 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 287 PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", number linear solve failures %" PetscInt_FMT " greater than current SNES allowed, stopping solve\n", snes->iter, snes->numLinearSolveFailures)); 288 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 289 break; 290 } 291 } 292 PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 293 snes->linear_its += lits; 294 PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 295 /* 296 if (snes->ops->precheck) { 297 PetscBool changed_y = PETSC_FALSE; 298 PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y); 299 } 300 301 if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W)); 302 */ 303 /* Compute a (scaled) negative update in the line search routine: 304 Y <- X - lambda*Y 305 and evaluate G = function(Y) (depends on the line search). 306 */ 307 PetscCall(VecCopy(Y, snes->vec_sol_update)); 308 ynorm = 1; 309 gnorm = vi->phinorm; 310 PetscCall(SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y)); 311 PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 312 PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm)); 313 PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)vi->phinorm, (double)gnorm, (double)ynorm, (int)lssucceed)); 314 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 315 if (snes->domainerror) { 316 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 317 sdm->ops->computefunction = vi->computeuserfunction; 318 PetscFunctionReturn(0); 319 } 320 if (lssucceed) { 321 if (++snes->numFailures >= snes->maxFailures) { 322 PetscBool ismin; 323 snes->reason = SNES_DIVERGED_LINE_SEARCH; 324 PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, vi->phi, X, gnorm, &ismin)); 325 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 326 break; 327 } 328 } 329 /* Update function and solution vectors */ 330 vi->phinorm = gnorm; 331 vi->merit = 0.5 * vi->phinorm * vi->phinorm; 332 /* Monitor convergence */ 333 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 334 snes->iter = i + 1; 335 snes->norm = vi->phinorm; 336 snes->xnorm = xnorm; 337 snes->ynorm = ynorm; 338 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 339 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 340 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 341 /* Test for convergence, xnorm = || X || */ 342 if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 343 PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, vi->phinorm, &snes->reason, snes->cnvP); 344 if (snes->reason) break; 345 } 346 if (i == maxits) { 347 PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 348 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 349 } 350 sdm->ops->computefunction = vi->computeuserfunction; 351 PetscFunctionReturn(0); 352 } 353 354 /* 355 SNESSetUp_VINEWTONSSLS - Sets up the internal data structures for the later use 356 of the SNES nonlinear solver. 357 358 Input Parameter: 359 . snes - the SNES context 360 361 Application Interface Routine: SNESSetUp() 362 363 Note: 364 For basic use of the SNES solvers, the user need not explicitly call 365 SNESSetUp(), since these actions will automatically occur during 366 the call to SNESSolve(). 367 */ 368 PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes) { 369 SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 370 371 PetscFunctionBegin; 372 PetscCall(SNESSetUp_VI(snes)); 373 PetscCall(VecDuplicate(snes->vec_sol, &vi->dpsi)); 374 PetscCall(VecDuplicate(snes->vec_sol, &vi->phi)); 375 PetscCall(VecDuplicate(snes->vec_sol, &vi->Da)); 376 PetscCall(VecDuplicate(snes->vec_sol, &vi->Db)); 377 PetscCall(VecDuplicate(snes->vec_sol, &vi->z)); 378 PetscCall(VecDuplicate(snes->vec_sol, &vi->t)); 379 PetscFunctionReturn(0); 380 } 381 382 PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes) { 383 SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 384 385 PetscFunctionBegin; 386 PetscCall(SNESReset_VI(snes)); 387 PetscCall(VecDestroy(&vi->dpsi)); 388 PetscCall(VecDestroy(&vi->phi)); 389 PetscCall(VecDestroy(&vi->Da)); 390 PetscCall(VecDestroy(&vi->Db)); 391 PetscCall(VecDestroy(&vi->z)); 392 PetscCall(VecDestroy(&vi->t)); 393 PetscFunctionReturn(0); 394 } 395 396 /* 397 SNESSetFromOptions_VINEWTONSSLS - Sets various parameters for the SNESVI method. 398 399 Input Parameter: 400 . snes - the SNES context 401 402 Application Interface Routine: SNESSetFromOptions() 403 */ 404 static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(SNES snes, PetscOptionItems *PetscOptionsObject) { 405 PetscFunctionBegin; 406 PetscCall(SNESSetFromOptions_VI(snes, PetscOptionsObject)); 407 PetscOptionsHeadBegin(PetscOptionsObject, "SNES semismooth method options"); 408 PetscOptionsHeadEnd(); 409 PetscFunctionReturn(0); 410 } 411 412 /*MC 413 SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method 414 415 Options Database Keys: 416 + -snes_type <vinewtonssls,vinewtonrsls> a semi-smooth solver, a reduced space active set method 417 - -snes_vi_monitor - prints the number of active constraints at each iteration. 418 419 Level: beginner 420 421 References: 422 + * - T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth 423 algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001). 424 - * - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale 425 Applications, Optimization Methods and Software, 21 (2006). 426 427 Notes: 428 This family of algorithm is much like an interior point method. 429 430 The reduced space active set solvers `SNESVINEWTONRSLS` provide an alternative approach that does not result in extremely ill-conditioned linear systems 431 432 .seealso: `SNESVINEWTONRSLS`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONRSLS`, `SNESNEWTONTR`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` 433 M*/ 434 PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes) { 435 SNES_VINEWTONSSLS *vi; 436 SNESLineSearch linesearch; 437 438 PetscFunctionBegin; 439 snes->ops->reset = SNESReset_VINEWTONSSLS; 440 snes->ops->setup = SNESSetUp_VINEWTONSSLS; 441 snes->ops->solve = SNESSolve_VINEWTONSSLS; 442 snes->ops->destroy = SNESDestroy_VI; 443 snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS; 444 snes->ops->view = NULL; 445 446 snes->usesksp = PETSC_TRUE; 447 snes->usesnpc = PETSC_FALSE; 448 449 PetscCall(SNESGetLineSearch(snes, &linesearch)); 450 if (!((PetscObject)linesearch)->type_name) { 451 PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 452 PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0)); 453 } 454 455 snes->alwayscomputesfinalresidual = PETSC_FALSE; 456 457 PetscCall(PetscNew(&vi)); 458 snes->data = (void *)vi; 459 460 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI)); 461 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI)); 462 PetscFunctionReturn(0); 463 } 464