1 2 #include <../src/snes/impls/vi/rs/virsimpl.h> /*I "petscsnes.h" I*/ 3 #include <petsc/private/dmimpl.h> 4 #include <petsc/private/vecimpl.h> 5 6 /*@ 7 SNESVIGetInactiveSet - Gets the global indices for the inactive set variables (these correspond to the degrees of freedom the linear 8 system is solved on) 9 10 Input Parameter: 11 . snes - the `SNES` context 12 13 Output Parameter: 14 . inact - inactive set index set 15 16 Level: advanced 17 18 .seealso: `SNESVINEWTONRSLS` 19 @*/ 20 PetscErrorCode SNESVIGetInactiveSet(SNES snes, IS *inact) 21 { 22 SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data; 23 24 PetscFunctionBegin; 25 *inact = vi->IS_inact; 26 PetscFunctionReturn(PETSC_SUCCESS); 27 } 28 29 /* 30 Provides a wrapper to a DM to allow it to be used to generated the interpolation/restriction from the DM for the smaller matrices and vectors 31 defined by the reduced space method. 32 33 Simple calls the regular DM interpolation and restricts it to operation on the variables not associated with active constraints. 34 */ 35 typedef struct { 36 PetscInt n; /* size of vectors in the reduced DM space */ 37 IS inactive; 38 39 PetscErrorCode (*createinterpolation)(DM, DM, Mat *, Vec *); /* DM's original routines */ 40 PetscErrorCode (*coarsen)(DM, MPI_Comm, DM *); 41 PetscErrorCode (*createglobalvector)(DM, Vec *); 42 PetscErrorCode (*createinjection)(DM, DM, Mat *); 43 PetscErrorCode (*hascreateinjection)(DM, PetscBool *); 44 45 DM dm; /* when destroying this object we need to reset the above function into the base DM */ 46 } DM_SNESVI; 47 48 /* 49 DMCreateGlobalVector_SNESVI - Creates global vector of the size of the reduced space 50 */ 51 static PetscErrorCode DMCreateGlobalVector_SNESVI(DM dm, Vec *vec) 52 { 53 PetscContainer isnes; 54 DM_SNESVI *dmsnesvi; 55 56 PetscFunctionBegin; 57 PetscCall(PetscObjectQuery((PetscObject)dm, "VI", (PetscObject *)&isnes)); 58 PetscCheck(isnes, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Composed SNES is missing"); 59 PetscCall(PetscContainerGetPointer(isnes, (void **)&dmsnesvi)); 60 PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)dm), dmsnesvi->n, PETSC_DETERMINE, vec)); 61 PetscCall(VecSetDM(*vec, dm)); 62 PetscFunctionReturn(PETSC_SUCCESS); 63 } 64 65 /* 66 DMCreateInterpolation_SNESVI - Modifieds the interpolation obtained from the DM by removing all rows and columns associated with active constraints. 67 */ 68 static PetscErrorCode DMCreateInterpolation_SNESVI(DM dm1, DM dm2, Mat *mat, Vec *vec) 69 { 70 PetscContainer isnes; 71 DM_SNESVI *dmsnesvi1, *dmsnesvi2; 72 Mat interp; 73 74 PetscFunctionBegin; 75 PetscCall(PetscObjectQuery((PetscObject)dm1, "VI", (PetscObject *)&isnes)); 76 PetscCheck(isnes, PetscObjectComm((PetscObject)dm1), PETSC_ERR_PLIB, "Composed VI data structure is missing"); 77 PetscCall(PetscContainerGetPointer(isnes, (void **)&dmsnesvi1)); 78 PetscCall(PetscObjectQuery((PetscObject)dm2, "VI", (PetscObject *)&isnes)); 79 PetscCheck(isnes, PetscObjectComm((PetscObject)dm2), PETSC_ERR_PLIB, "Composed VI data structure is missing"); 80 PetscCall(PetscContainerGetPointer(isnes, (void **)&dmsnesvi2)); 81 82 PetscCall((*dmsnesvi1->createinterpolation)(dm1, dm2, &interp, NULL)); 83 PetscCall(MatCreateSubMatrix(interp, dmsnesvi2->inactive, dmsnesvi1->inactive, MAT_INITIAL_MATRIX, mat)); 84 PetscCall(MatDestroy(&interp)); 85 *vec = NULL; 86 PetscFunctionReturn(PETSC_SUCCESS); 87 } 88 89 /* 90 DMCoarsen_SNESVI - Computes the regular coarsened DM then computes additional information about its inactive set 91 */ 92 static PetscErrorCode DMCoarsen_SNESVI(DM dm1, MPI_Comm comm, DM *dm2) 93 { 94 PetscContainer isnes; 95 DM_SNESVI *dmsnesvi1; 96 Vec finemarked, coarsemarked; 97 IS inactive; 98 Mat inject; 99 const PetscInt *index; 100 PetscInt n, k, cnt = 0, rstart, *coarseindex; 101 PetscScalar *marked; 102 103 PetscFunctionBegin; 104 PetscCall(PetscObjectQuery((PetscObject)dm1, "VI", (PetscObject *)&isnes)); 105 PetscCheck(isnes, PetscObjectComm((PetscObject)dm1), PETSC_ERR_PLIB, "Composed VI data structure is missing"); 106 PetscCall(PetscContainerGetPointer(isnes, (void **)&dmsnesvi1)); 107 108 /* get the original coarsen */ 109 PetscCall((*dmsnesvi1->coarsen)(dm1, comm, dm2)); 110 111 /* not sure why this extra reference is needed, but without the dm2 disappears too early */ 112 /* Updating the KSPCreateVecs() to avoid using DMGetGlobalVector() when matrix is available removes the need for this reference? */ 113 /* PetscCall(PetscObjectReference((PetscObject)*dm2));*/ 114 115 /* need to set back global vectors in order to use the original injection */ 116 PetscCall(DMClearGlobalVectors(dm1)); 117 118 dm1->ops->createglobalvector = dmsnesvi1->createglobalvector; 119 120 PetscCall(DMCreateGlobalVector(dm1, &finemarked)); 121 PetscCall(DMCreateGlobalVector(*dm2, &coarsemarked)); 122 123 /* 124 fill finemarked with locations of inactive points 125 */ 126 PetscCall(ISGetIndices(dmsnesvi1->inactive, &index)); 127 PetscCall(ISGetLocalSize(dmsnesvi1->inactive, &n)); 128 PetscCall(VecSet(finemarked, 0.0)); 129 for (k = 0; k < n; k++) PetscCall(VecSetValue(finemarked, index[k], 1.0, INSERT_VALUES)); 130 PetscCall(VecAssemblyBegin(finemarked)); 131 PetscCall(VecAssemblyEnd(finemarked)); 132 133 PetscCall(DMCreateInjection(*dm2, dm1, &inject)); 134 PetscCall(MatRestrict(inject, finemarked, coarsemarked)); 135 PetscCall(MatDestroy(&inject)); 136 137 /* 138 create index set list of coarse inactive points from coarsemarked 139 */ 140 PetscCall(VecGetLocalSize(coarsemarked, &n)); 141 PetscCall(VecGetOwnershipRange(coarsemarked, &rstart, NULL)); 142 PetscCall(VecGetArray(coarsemarked, &marked)); 143 for (k = 0; k < n; k++) { 144 if (marked[k] != 0.0) cnt++; 145 } 146 PetscCall(PetscMalloc1(cnt, &coarseindex)); 147 cnt = 0; 148 for (k = 0; k < n; k++) { 149 if (marked[k] != 0.0) coarseindex[cnt++] = k + rstart; 150 } 151 PetscCall(VecRestoreArray(coarsemarked, &marked)); 152 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)coarsemarked), cnt, coarseindex, PETSC_OWN_POINTER, &inactive)); 153 154 PetscCall(DMClearGlobalVectors(dm1)); 155 156 dm1->ops->createglobalvector = DMCreateGlobalVector_SNESVI; 157 158 PetscCall(DMSetVI(*dm2, inactive)); 159 160 PetscCall(VecDestroy(&finemarked)); 161 PetscCall(VecDestroy(&coarsemarked)); 162 PetscCall(ISDestroy(&inactive)); 163 PetscFunctionReturn(PETSC_SUCCESS); 164 } 165 166 static PetscErrorCode DMDestroy_SNESVI(DM_SNESVI *dmsnesvi) 167 { 168 PetscFunctionBegin; 169 /* reset the base methods in the DM object that were changed when the DM_SNESVI was reset */ 170 dmsnesvi->dm->ops->createinterpolation = dmsnesvi->createinterpolation; 171 dmsnesvi->dm->ops->coarsen = dmsnesvi->coarsen; 172 dmsnesvi->dm->ops->createglobalvector = dmsnesvi->createglobalvector; 173 dmsnesvi->dm->ops->createinjection = dmsnesvi->createinjection; 174 dmsnesvi->dm->ops->hascreateinjection = dmsnesvi->hascreateinjection; 175 /* need to clear out this vectors because some of them may not have a reference to the DM 176 but they are counted as having references to the DM in DMDestroy() */ 177 PetscCall(DMClearGlobalVectors(dmsnesvi->dm)); 178 179 PetscCall(ISDestroy(&dmsnesvi->inactive)); 180 PetscCall(PetscFree(dmsnesvi)); 181 PetscFunctionReturn(PETSC_SUCCESS); 182 } 183 184 /*@ 185 DMSetVI - Marks a DM as associated with a VI problem. This causes the interpolation/restriction operators to 186 be restricted to only those variables NOT associated with active constraints. 187 188 Logically Collective 189 190 Input Parameters: 191 + dm - the `DM` object 192 - inactive - an `IS` indicating which points are currently not active 193 194 Level: intermediate 195 196 .seealso: `SNESVINEWTONRSLS`, `SNESVIGetInactiveSet()` 197 @*/ 198 PetscErrorCode DMSetVI(DM dm, IS inactive) 199 { 200 PetscContainer isnes; 201 DM_SNESVI *dmsnesvi; 202 203 PetscFunctionBegin; 204 if (!dm) PetscFunctionReturn(PETSC_SUCCESS); 205 206 PetscCall(PetscObjectReference((PetscObject)inactive)); 207 208 PetscCall(PetscObjectQuery((PetscObject)dm, "VI", (PetscObject *)&isnes)); 209 if (!isnes) { 210 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)dm), &isnes)); 211 PetscCall(PetscContainerSetUserDestroy(isnes, (PetscErrorCode(*)(void *))DMDestroy_SNESVI)); 212 PetscCall(PetscNew(&dmsnesvi)); 213 PetscCall(PetscContainerSetPointer(isnes, (void *)dmsnesvi)); 214 PetscCall(PetscObjectCompose((PetscObject)dm, "VI", (PetscObject)isnes)); 215 PetscCall(PetscContainerDestroy(&isnes)); 216 217 dmsnesvi->createinterpolation = dm->ops->createinterpolation; 218 dm->ops->createinterpolation = DMCreateInterpolation_SNESVI; 219 dmsnesvi->coarsen = dm->ops->coarsen; 220 dm->ops->coarsen = DMCoarsen_SNESVI; 221 dmsnesvi->createglobalvector = dm->ops->createglobalvector; 222 dm->ops->createglobalvector = DMCreateGlobalVector_SNESVI; 223 dmsnesvi->createinjection = dm->ops->createinjection; 224 dm->ops->createinjection = NULL; 225 dmsnesvi->hascreateinjection = dm->ops->hascreateinjection; 226 dm->ops->hascreateinjection = NULL; 227 } else { 228 PetscCall(PetscContainerGetPointer(isnes, (void **)&dmsnesvi)); 229 PetscCall(ISDestroy(&dmsnesvi->inactive)); 230 } 231 PetscCall(DMClearGlobalVectors(dm)); 232 PetscCall(ISGetLocalSize(inactive, &dmsnesvi->n)); 233 234 dmsnesvi->inactive = inactive; 235 dmsnesvi->dm = dm; 236 PetscFunctionReturn(PETSC_SUCCESS); 237 } 238 239 /* 240 DMDestroyVI - Frees the DM_SNESVI object contained in the DM 241 - also resets the function pointers in the DM for createinterpolation() etc to use the original DM 242 */ 243 PetscErrorCode DMDestroyVI(DM dm) 244 { 245 PetscFunctionBegin; 246 if (!dm) PetscFunctionReturn(PETSC_SUCCESS); 247 PetscCall(PetscObjectCompose((PetscObject)dm, "VI", (PetscObject)NULL)); 248 PetscFunctionReturn(PETSC_SUCCESS); 249 } 250 251 /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */ 252 static PetscErrorCode SNESCreateSubVectors_VINEWTONRSLS(SNES snes, PetscInt n, Vec *newv) 253 { 254 Vec v; 255 256 PetscFunctionBegin; 257 PetscCall(VecCreate(PetscObjectComm((PetscObject)snes), &v)); 258 PetscCall(VecSetSizes(v, n, PETSC_DECIDE)); 259 PetscCall(VecSetType(v, VECSTANDARD)); 260 *newv = v; 261 PetscFunctionReturn(PETSC_SUCCESS); 262 } 263 264 /* Resets the snes PC and KSP when the active set sizes change */ 265 static PetscErrorCode SNESVIResetPCandKSP(SNES snes, Mat Amat, Mat Pmat) 266 { 267 KSP snesksp; 268 269 PetscFunctionBegin; 270 PetscCall(SNESGetKSP(snes, &snesksp)); 271 PetscCall(KSPReset(snesksp)); 272 PetscCall(KSPResetFromOptions(snesksp)); 273 274 /* 275 KSP kspnew; 276 PC pcnew; 277 MatSolverType stype; 278 279 PetscCall(KSPCreate(PetscObjectComm((PetscObject)snes),&kspnew)); 280 kspnew->pc_side = snesksp->pc_side; 281 kspnew->rtol = snesksp->rtol; 282 kspnew->abstol = snesksp->abstol; 283 kspnew->max_it = snesksp->max_it; 284 PetscCall(KSPSetType(kspnew,((PetscObject)snesksp)->type_name)); 285 PetscCall(KSPGetPC(kspnew,&pcnew)); 286 PetscCall(PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name)); 287 PetscCall(PCSetOperators(kspnew->pc,Amat,Pmat)); 288 PetscCall(PCFactorGetMatSolverType(snesksp->pc,&stype)); 289 PetscCall(PCFactorSetMatSolverType(kspnew->pc,stype)); 290 PetscCall(KSPDestroy(&snesksp)); 291 snes->ksp = kspnew; 292 PetscCall(KSPSetFromOptions(kspnew));*/ 293 PetscFunctionReturn(PETSC_SUCCESS); 294 } 295 296 /* Variational Inequality solver using reduce space method. No semismooth algorithm is 297 implemented in this algorithm. It basically identifies the active constraints and does 298 a linear solve on the other variables (those not associated with the active constraints). */ 299 static PetscErrorCode SNESSolve_VINEWTONRSLS(SNES snes) 300 { 301 SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data; 302 PetscInt maxits, i, lits; 303 SNESLineSearchReason lssucceed; 304 PetscReal fnorm, gnorm, xnorm = 0, ynorm; 305 Vec Y, X, F; 306 KSPConvergedReason kspreason; 307 KSP ksp; 308 PC pc; 309 310 PetscFunctionBegin; 311 /* Multigrid must use Galerkin for coarse grids with active set/reduced space methods; cannot rediscretize on coarser grids*/ 312 PetscCall(SNESGetKSP(snes, &ksp)); 313 PetscCall(KSPGetPC(ksp, &pc)); 314 PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_BOTH)); 315 316 snes->numFailures = 0; 317 snes->numLinearSolveFailures = 0; 318 snes->reason = SNES_CONVERGED_ITERATING; 319 320 maxits = snes->max_its; /* maximum number of iterations */ 321 X = snes->vec_sol; /* solution vector */ 322 F = snes->vec_func; /* residual vector */ 323 Y = snes->work[0]; /* work vectors */ 324 325 PetscCall(SNESLineSearchSetVIFunctions(snes->linesearch, SNESVIProjectOntoBounds, SNESVIComputeInactiveSetFnorm)); 326 PetscCall(SNESLineSearchSetVecs(snes->linesearch, X, NULL, NULL, NULL, NULL)); 327 PetscCall(SNESLineSearchSetUp(snes->linesearch)); 328 329 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 330 snes->iter = 0; 331 snes->norm = 0.0; 332 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 333 334 PetscCall(SNESVIProjectOntoBounds(snes, X)); 335 PetscCall(SNESComputeFunction(snes, X, F)); 336 PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm)); 337 PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- ||x|| */ 338 SNESCheckFunctionNorm(snes, fnorm); 339 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 340 snes->norm = fnorm; 341 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 342 PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 343 344 /* test convergence */ 345 PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 346 PetscCall(SNESMonitor(snes, 0, fnorm)); 347 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 348 349 for (i = 0; i < maxits; i++) { 350 IS IS_act; /* _act -> active set _inact -> inactive set */ 351 IS IS_redact; /* redundant active set */ 352 VecScatter scat_act, scat_inact; 353 PetscInt nis_act, nis_inact; 354 Vec Y_act, Y_inact, F_inact; 355 Mat jac_inact_inact, prejac_inact_inact; 356 PetscBool isequal; 357 358 /* Call general purpose update function */ 359 PetscTryTypeMethod(snes, update, snes->iter); 360 PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 361 SNESCheckJacobianDomainerror(snes); 362 363 /* Create active and inactive index sets */ 364 365 /*original 366 PetscCall(SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&vi->IS_inact)); 367 */ 368 PetscCall(SNESVIGetActiveSetIS(snes, X, F, &IS_act)); 369 370 if (vi->checkredundancy) { 371 PetscCall((*vi->checkredundancy)(snes, IS_act, &IS_redact, vi->ctxP)); 372 if (IS_redact) { 373 PetscCall(ISSort(IS_redact)); 374 PetscCall(ISComplement(IS_redact, X->map->rstart, X->map->rend, &vi->IS_inact)); 375 PetscCall(ISDestroy(&IS_redact)); 376 } else { 377 PetscCall(ISComplement(IS_act, X->map->rstart, X->map->rend, &vi->IS_inact)); 378 } 379 } else { 380 PetscCall(ISComplement(IS_act, X->map->rstart, X->map->rend, &vi->IS_inact)); 381 } 382 383 /* Create inactive set submatrix */ 384 PetscCall(MatCreateSubMatrix(snes->jacobian, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &jac_inact_inact)); 385 386 if (0) { /* Dead code (temporary developer hack) */ 387 IS keptrows; 388 PetscCall(MatFindNonzeroRows(jac_inact_inact, &keptrows)); 389 if (keptrows) { 390 PetscInt cnt, *nrows, k; 391 const PetscInt *krows, *inact; 392 PetscInt rstart; 393 394 PetscCall(MatGetOwnershipRange(jac_inact_inact, &rstart, NULL)); 395 PetscCall(MatDestroy(&jac_inact_inact)); 396 PetscCall(ISDestroy(&IS_act)); 397 398 PetscCall(ISGetLocalSize(keptrows, &cnt)); 399 PetscCall(ISGetIndices(keptrows, &krows)); 400 PetscCall(ISGetIndices(vi->IS_inact, &inact)); 401 PetscCall(PetscMalloc1(cnt, &nrows)); 402 for (k = 0; k < cnt; k++) nrows[k] = inact[krows[k] - rstart]; 403 PetscCall(ISRestoreIndices(keptrows, &krows)); 404 PetscCall(ISRestoreIndices(vi->IS_inact, &inact)); 405 PetscCall(ISDestroy(&keptrows)); 406 PetscCall(ISDestroy(&vi->IS_inact)); 407 408 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), cnt, nrows, PETSC_OWN_POINTER, &vi->IS_inact)); 409 PetscCall(ISComplement(vi->IS_inact, F->map->rstart, F->map->rend, &IS_act)); 410 PetscCall(MatCreateSubMatrix(snes->jacobian, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &jac_inact_inact)); 411 } 412 } 413 PetscCall(DMSetVI(snes->dm, vi->IS_inact)); 414 /* remove later */ 415 416 /* 417 PetscCall(VecView(vi->xu,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xu))->comm))); 418 PetscCall(VecView(vi->xl,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xl))->comm))); 419 PetscCall(VecView(X,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)X)))); 420 PetscCall(VecView(F,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)F)))); 421 PetscCall(ISView(vi->IS_inact,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)vi->IS_inact)))); 422 */ 423 424 /* Get sizes of active and inactive sets */ 425 PetscCall(ISGetLocalSize(IS_act, &nis_act)); 426 PetscCall(ISGetLocalSize(vi->IS_inact, &nis_inact)); 427 428 /* Create active and inactive set vectors */ 429 PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_inact, &F_inact)); 430 PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_act, &Y_act)); 431 PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_inact, &Y_inact)); 432 433 /* Create scatter contexts */ 434 PetscCall(VecScatterCreate(Y, IS_act, Y_act, NULL, &scat_act)); 435 PetscCall(VecScatterCreate(Y, vi->IS_inact, Y_inact, NULL, &scat_inact)); 436 437 /* Do a vec scatter to active and inactive set vectors */ 438 PetscCall(VecScatterBegin(scat_inact, F, F_inact, INSERT_VALUES, SCATTER_FORWARD)); 439 PetscCall(VecScatterEnd(scat_inact, F, F_inact, INSERT_VALUES, SCATTER_FORWARD)); 440 441 PetscCall(VecScatterBegin(scat_act, Y, Y_act, INSERT_VALUES, SCATTER_FORWARD)); 442 PetscCall(VecScatterEnd(scat_act, Y, Y_act, INSERT_VALUES, SCATTER_FORWARD)); 443 444 PetscCall(VecScatterBegin(scat_inact, Y, Y_inact, INSERT_VALUES, SCATTER_FORWARD)); 445 PetscCall(VecScatterEnd(scat_inact, Y, Y_inact, INSERT_VALUES, SCATTER_FORWARD)); 446 447 /* Active set direction = 0 */ 448 PetscCall(VecSet(Y_act, 0)); 449 if (snes->jacobian != snes->jacobian_pre) { 450 PetscCall(MatCreateSubMatrix(snes->jacobian_pre, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &prejac_inact_inact)); 451 } else prejac_inact_inact = jac_inact_inact; 452 453 PetscCall(ISEqual(vi->IS_inact_prev, vi->IS_inact, &isequal)); 454 if (!isequal) { 455 PetscCall(SNESVIResetPCandKSP(snes, jac_inact_inact, prejac_inact_inact)); 456 PetscCall(PCFieldSplitRestrictIS(pc, vi->IS_inact)); 457 } 458 459 /* PetscCall(ISView(vi->IS_inact,0)); */ 460 /* PetscCall(ISView(IS_act,0));*/ 461 /* ierr = MatView(snes->jacobian_pre,0); */ 462 463 PetscCall(KSPSetOperators(snes->ksp, jac_inact_inact, prejac_inact_inact)); 464 PetscCall(KSPSetUp(snes->ksp)); 465 { 466 PC pc; 467 PetscBool flg; 468 PetscCall(KSPGetPC(snes->ksp, &pc)); 469 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &flg)); 470 if (flg) { 471 KSP *subksps; 472 PetscCall(PCFieldSplitGetSubKSP(pc, NULL, &subksps)); 473 PetscCall(KSPGetPC(subksps[0], &pc)); 474 PetscCall(PetscFree(subksps)); 475 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &flg)); 476 if (flg) { 477 PetscInt n, N = 101 * 101, j, cnts[3] = {0, 0, 0}; 478 const PetscInt *ii; 479 480 PetscCall(ISGetSize(vi->IS_inact, &n)); 481 PetscCall(ISGetIndices(vi->IS_inact, &ii)); 482 for (j = 0; j < n; j++) { 483 if (ii[j] < N) cnts[0]++; 484 else if (ii[j] < 2 * N) cnts[1]++; 485 else if (ii[j] < 3 * N) cnts[2]++; 486 } 487 PetscCall(ISRestoreIndices(vi->IS_inact, &ii)); 488 489 PetscCall(PCBJacobiSetTotalBlocks(pc, 3, cnts)); 490 } 491 } 492 } 493 494 PetscCall(KSPSolve(snes->ksp, F_inact, Y_inact)); 495 PetscCall(VecScatterBegin(scat_act, Y_act, Y, INSERT_VALUES, SCATTER_REVERSE)); 496 PetscCall(VecScatterEnd(scat_act, Y_act, Y, INSERT_VALUES, SCATTER_REVERSE)); 497 PetscCall(VecScatterBegin(scat_inact, Y_inact, Y, INSERT_VALUES, SCATTER_REVERSE)); 498 PetscCall(VecScatterEnd(scat_inact, Y_inact, Y, INSERT_VALUES, SCATTER_REVERSE)); 499 500 PetscCall(VecDestroy(&F_inact)); 501 PetscCall(VecDestroy(&Y_act)); 502 PetscCall(VecDestroy(&Y_inact)); 503 PetscCall(VecScatterDestroy(&scat_act)); 504 PetscCall(VecScatterDestroy(&scat_inact)); 505 PetscCall(ISDestroy(&IS_act)); 506 if (!isequal) { 507 PetscCall(ISDestroy(&vi->IS_inact_prev)); 508 PetscCall(ISDuplicate(vi->IS_inact, &vi->IS_inact_prev)); 509 } 510 PetscCall(ISDestroy(&vi->IS_inact)); 511 PetscCall(MatDestroy(&jac_inact_inact)); 512 if (snes->jacobian != snes->jacobian_pre) PetscCall(MatDestroy(&prejac_inact_inact)); 513 514 PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason)); 515 if (kspreason < 0) { 516 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 517 PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", number linear solve failures %" PetscInt_FMT " greater than current SNES allowed, stopping solve\n", snes->iter, snes->numLinearSolveFailures)); 518 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 519 break; 520 } 521 } 522 523 PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 524 snes->linear_its += lits; 525 PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 526 /* 527 if (snes->ops->precheck) { 528 PetscBool changed_y = PETSC_FALSE; 529 PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y); 530 } 531 532 if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W)); 533 */ 534 /* Compute a (scaled) negative update in the line search routine: 535 Y <- X - lambda*Y 536 and evaluate G = function(Y) (depends on the line search). 537 */ 538 PetscCall(VecCopy(Y, snes->vec_sol_update)); 539 ynorm = 1; 540 gnorm = fnorm; 541 PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &gnorm, Y)); 542 PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 543 PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm)); 544 PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lssucceed)); 545 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 546 if (snes->domainerror) { 547 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 548 PetscCall(DMDestroyVI(snes->dm)); 549 PetscFunctionReturn(PETSC_SUCCESS); 550 } 551 if (lssucceed) { 552 if (++snes->numFailures >= snes->maxFailures) { 553 PetscBool ismin; 554 snes->reason = SNES_DIVERGED_LINE_SEARCH; 555 PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, F, X, gnorm, &ismin)); 556 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 557 break; 558 } 559 } 560 PetscCall(DMDestroyVI(snes->dm)); 561 /* Update function and solution vectors */ 562 fnorm = gnorm; 563 /* Monitor convergence */ 564 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 565 snes->iter = i + 1; 566 snes->norm = fnorm; 567 snes->xnorm = xnorm; 568 snes->ynorm = ynorm; 569 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 570 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 571 /* Test for convergence, xnorm = || X || */ 572 if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 573 PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 574 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 575 if (snes->reason) break; 576 } 577 /* make sure that the VI information attached to the DM is removed if the for loop above was broken early due to some exceptional conditional */ 578 PetscCall(DMDestroyVI(snes->dm)); 579 PetscFunctionReturn(PETSC_SUCCESS); 580 } 581 582 /*@C 583 SNESVISetRedundancyCheck - Provide a function to check for any redundancy in the VI active set 584 585 Logically Collective 586 587 Input Parameters: 588 + snes - the `SNESVINEWTONRSLS` context 589 . func - the function to check of redundancies 590 - ctx - optional context used by the function 591 592 Level: advanced 593 594 Note: 595 Sometimes the inactive set will result in a non-singular sub-Jacobian problem that needs to be solved, this allows the user, 596 when they know more about their specific problem to provide a function that removes the redundancy that results in the singular linear system 597 598 .seealso: `SNESVINEWTONRSLS`, `SNESVIGetInactiveSet()`, `DMSetVI()` 599 @*/ 600 PetscErrorCode SNESVISetRedundancyCheck(SNES snes, PetscErrorCode (*func)(SNES, IS, IS *, void *), void *ctx) 601 { 602 SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data; 603 604 PetscFunctionBegin; 605 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 606 vi->checkredundancy = func; 607 vi->ctxP = ctx; 608 PetscFunctionReturn(PETSC_SUCCESS); 609 } 610 611 #if defined(PETSC_HAVE_MATLAB) 612 #include <engine.h> 613 #include <mex.h> 614 typedef struct { 615 char *funcname; 616 mxArray *ctx; 617 } SNESMatlabContext; 618 619 PetscErrorCode SNESVIRedundancyCheck_Matlab(SNES snes, IS is_act, IS *is_redact, void *ctx) 620 { 621 SNESMatlabContext *sctx = (SNESMatlabContext *)ctx; 622 int nlhs = 1, nrhs = 5; 623 mxArray *plhs[1], *prhs[5]; 624 long long int l1 = 0, l2 = 0, ls = 0; 625 PetscInt *indices = NULL; 626 627 PetscFunctionBegin; 628 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 629 PetscValidHeaderSpecific(is_act, IS_CLASSID, 2); 630 PetscAssertPointer(is_redact, 3); 631 PetscCheckSameComm(snes, 1, is_act, 2); 632 633 /* Create IS for reduced active set of size 0, its size and indices will 634 bet set by the Matlab function */ 635 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), 0, indices, PETSC_OWN_POINTER, is_redact)); 636 /* call Matlab function in ctx */ 637 PetscCall(PetscArraycpy(&ls, &snes, 1)); 638 PetscCall(PetscArraycpy(&l1, &is_act, 1)); 639 PetscCall(PetscArraycpy(&l2, is_redact, 1)); 640 prhs[0] = mxCreateDoubleScalar((double)ls); 641 prhs[1] = mxCreateDoubleScalar((double)l1); 642 prhs[2] = mxCreateDoubleScalar((double)l2); 643 prhs[3] = mxCreateString(sctx->funcname); 644 prhs[4] = sctx->ctx; 645 PetscCall(mexCallMATLAB(nlhs, plhs, nrhs, prhs, "PetscSNESVIRedundancyCheckInternal")); 646 PetscCall(mxGetScalar(plhs[0])); 647 mxDestroyArray(prhs[0]); 648 mxDestroyArray(prhs[1]); 649 mxDestroyArray(prhs[2]); 650 mxDestroyArray(prhs[3]); 651 mxDestroyArray(plhs[0]); 652 PetscFunctionReturn(PETSC_SUCCESS); 653 } 654 655 PetscErrorCode SNESVISetRedundancyCheckMatlab(SNES snes, const char *func, mxArray *ctx) 656 { 657 SNESMatlabContext *sctx; 658 659 PetscFunctionBegin; 660 /* currently sctx is memory bleed */ 661 PetscCall(PetscNew(&sctx)); 662 PetscCall(PetscStrallocpy(func, &sctx->funcname)); 663 sctx->ctx = mxDuplicateArray(ctx); 664 PetscCall(SNESVISetRedundancyCheck(snes, SNESVIRedundancyCheck_Matlab, sctx)); 665 PetscFunctionReturn(PETSC_SUCCESS); 666 } 667 668 #endif 669 670 static PetscErrorCode SNESSetUp_VINEWTONRSLS(SNES snes) 671 { 672 SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data; 673 PetscInt *indices; 674 PetscInt i, n, rstart, rend; 675 SNESLineSearch linesearch; 676 677 PetscFunctionBegin; 678 PetscCall(SNESSetUp_VI(snes)); 679 680 /* Set up previous active index set for the first snes solve 681 vi->IS_inact_prev = 0,1,2,....N */ 682 683 PetscCall(VecGetOwnershipRange(snes->vec_sol, &rstart, &rend)); 684 PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 685 PetscCall(PetscMalloc1(n, &indices)); 686 for (i = 0; i < n; i++) indices[i] = rstart + i; 687 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), n, indices, PETSC_OWN_POINTER, &vi->IS_inact_prev)); 688 689 /* set the line search functions */ 690 if (!snes->linesearch) { 691 PetscCall(SNESGetLineSearch(snes, &linesearch)); 692 PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 693 } 694 PetscFunctionReturn(PETSC_SUCCESS); 695 } 696 697 static PetscErrorCode SNESReset_VINEWTONRSLS(SNES snes) 698 { 699 SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data; 700 701 PetscFunctionBegin; 702 PetscCall(SNESReset_VI(snes)); 703 PetscCall(ISDestroy(&vi->IS_inact_prev)); 704 PetscFunctionReturn(PETSC_SUCCESS); 705 } 706 707 /*MC 708 SNESVINEWTONRSLS - Reduced space active set solvers for variational inequalities based on Newton's method 709 710 Options Database Keys: 711 + -snes_type <vinewtonssls,vinewtonrsls> - a semi-smooth solver, a reduced space active set method 712 - -snes_vi_monitor - prints the number of active constraints at each iteration. 713 714 Level: beginner 715 716 References: 717 . * - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale 718 Applications, Optimization Methods and Software, 21 (2006). 719 720 Note: 721 At each set of this methods the algorithm produces an inactive set of variables that are constrained to their current values 722 (because changing these values would result in those variables no longer satisfying the inequality constraints) 723 and produces a step direction by solving the linear system arising from the Jacobian with the inactive variables removed. In other 724 words on a reduced space of the solution space. Based on the Newton update it then adjusts the inactive sep for the next iteration. 725 726 .seealso: `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONSSLS`, `SNESNEWTONTR`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESVIGetInactiveSet()`, `DMSetVI()`, `SNESVISetRedundancyCheck()` 727 M*/ 728 PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONRSLS(SNES snes) 729 { 730 SNES_VINEWTONRSLS *vi; 731 SNESLineSearch linesearch; 732 733 PetscFunctionBegin; 734 snes->ops->reset = SNESReset_VINEWTONRSLS; 735 snes->ops->setup = SNESSetUp_VINEWTONRSLS; 736 snes->ops->solve = SNESSolve_VINEWTONRSLS; 737 snes->ops->destroy = SNESDestroy_VI; 738 snes->ops->setfromoptions = SNESSetFromOptions_VI; 739 snes->ops->view = NULL; 740 snes->ops->converged = SNESConvergedDefault_VI; 741 742 snes->usesksp = PETSC_TRUE; 743 snes->usesnpc = PETSC_FALSE; 744 745 PetscCall(SNESGetLineSearch(snes, &linesearch)); 746 if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 747 PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0)); 748 749 snes->alwayscomputesfinalresidual = PETSC_TRUE; 750 751 PetscCall(PetscNew(&vi)); 752 snes->data = (void *)vi; 753 vi->checkredundancy = NULL; 754 755 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI)); 756 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI)); 757 PetscFunctionReturn(PETSC_SUCCESS); 758 } 759