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