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