1c2fc9fa9SBarry Smith 2c2fc9fa9SBarry Smith #include <../src/snes/impls/vi/rs/virsimpl.h> /*I "petscsnes.h" I*/ 3af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 4af0996ceSBarry Smith #include <petsc/private/vecimpl.h> 5c2fc9fa9SBarry Smith 6*f6dfbefdSBarry Smith /*@ 7c2fc9fa9SBarry Smith SNESVIGetInactiveSet - Gets the global indices for the inactive set variables (these correspond to the degrees of freedom the linear 8c2fc9fa9SBarry Smith system is solved on) 9c2fc9fa9SBarry Smith 107a7aea1fSJed Brown Input parameter: 11*f6dfbefdSBarry Smith . snes - the `SNES` context 12c2fc9fa9SBarry Smith 137a7aea1fSJed Brown Output parameter: 14d5f1b7e6SEd Bueler . inact - inactive set index set 15c2fc9fa9SBarry Smith 16*f6dfbefdSBarry Smith .seealso: `SNESVINEWTONRSLS` 17*f6dfbefdSBarry Smith @*/ 189371c9d4SSatish Balay PetscErrorCode SNESVIGetInactiveSet(SNES snes, IS *inact) { 19f450aa47SBarry Smith SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data; 206e111a19SKarl Rupp 21c2fc9fa9SBarry Smith PetscFunctionBegin; 22f009fc93SPatrick Farrell *inact = vi->IS_inact; 23c2fc9fa9SBarry Smith PetscFunctionReturn(0); 24c2fc9fa9SBarry Smith } 25c2fc9fa9SBarry Smith 26c2fc9fa9SBarry Smith /* 27c2fc9fa9SBarry Smith 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 28c2fc9fa9SBarry Smith defined by the reduced space method. 29c2fc9fa9SBarry Smith 30c2fc9fa9SBarry Smith Simple calls the regular DM interpolation and restricts it to operation on the variables not associated with active constraints. 31d49cd2b7SBarry Smith */ 32c2fc9fa9SBarry Smith typedef struct { 33c2fc9fa9SBarry Smith PetscInt n; /* size of vectors in the reduced DM space */ 34c2fc9fa9SBarry Smith IS inactive; 35f5af7f23SKarl Rupp 3625296bd5SBarry Smith PetscErrorCode (*createinterpolation)(DM, DM, Mat *, Vec *); /* DM's original routines */ 37c2fc9fa9SBarry Smith PetscErrorCode (*coarsen)(DM, MPI_Comm, DM *); 38c2fc9fa9SBarry Smith PetscErrorCode (*createglobalvector)(DM, Vec *); 395a84ad33SLisandro Dalcin PetscErrorCode (*createinjection)(DM, DM, Mat *); 404a7a4c06SLawrence Mitchell PetscErrorCode (*hascreateinjection)(DM, PetscBool *); 41f5af7f23SKarl Rupp 42c2fc9fa9SBarry Smith DM dm; /* when destroying this object we need to reset the above function into the base DM */ 43c2fc9fa9SBarry Smith } DM_SNESVI; 44c2fc9fa9SBarry Smith 45c2fc9fa9SBarry Smith /* 46c2fc9fa9SBarry Smith DMCreateGlobalVector_SNESVI - Creates global vector of the size of the reduced space 47c2fc9fa9SBarry Smith */ 489371c9d4SSatish Balay PetscErrorCode DMCreateGlobalVector_SNESVI(DM dm, Vec *vec) { 49c2fc9fa9SBarry Smith PetscContainer isnes; 50c2fc9fa9SBarry Smith DM_SNESVI *dmsnesvi; 51c2fc9fa9SBarry Smith 52c2fc9fa9SBarry Smith PetscFunctionBegin; 539566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "VI", (PetscObject *)&isnes)); 5428b400f6SJacob Faibussowitsch PetscCheck(isnes, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Composed SNES is missing"); 559566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(isnes, (void **)&dmsnesvi)); 569566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)dm), dmsnesvi->n, PETSC_DETERMINE, vec)); 579566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec, dm)); 58c2fc9fa9SBarry Smith PetscFunctionReturn(0); 59c2fc9fa9SBarry Smith } 60c2fc9fa9SBarry Smith 61c2fc9fa9SBarry Smith /* 62e727c939SJed Brown DMCreateInterpolation_SNESVI - Modifieds the interpolation obtained from the DM by removing all rows and columns associated with active constraints. 63c2fc9fa9SBarry Smith */ 649371c9d4SSatish Balay PetscErrorCode DMCreateInterpolation_SNESVI(DM dm1, DM dm2, Mat *mat, Vec *vec) { 65c2fc9fa9SBarry Smith PetscContainer isnes; 66c2fc9fa9SBarry Smith DM_SNESVI *dmsnesvi1, *dmsnesvi2; 67c2fc9fa9SBarry Smith Mat interp; 68c2fc9fa9SBarry Smith 69c2fc9fa9SBarry Smith PetscFunctionBegin; 709566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm1, "VI", (PetscObject *)&isnes)); 7128b400f6SJacob Faibussowitsch PetscCheck(isnes, PetscObjectComm((PetscObject)dm1), PETSC_ERR_PLIB, "Composed VI data structure is missing"); 729566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(isnes, (void **)&dmsnesvi1)); 739566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm2, "VI", (PetscObject *)&isnes)); 7428b400f6SJacob Faibussowitsch PetscCheck(isnes, PetscObjectComm((PetscObject)dm2), PETSC_ERR_PLIB, "Composed VI data structure is missing"); 759566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(isnes, (void **)&dmsnesvi2)); 76c2fc9fa9SBarry Smith 779566063dSJacob Faibussowitsch PetscCall((*dmsnesvi1->createinterpolation)(dm1, dm2, &interp, NULL)); 789566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(interp, dmsnesvi2->inactive, dmsnesvi1->inactive, MAT_INITIAL_MATRIX, mat)); 799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&interp)); 809e5d0892SLisandro Dalcin *vec = NULL; 81c2fc9fa9SBarry Smith PetscFunctionReturn(0); 82c2fc9fa9SBarry Smith } 83c2fc9fa9SBarry Smith 84c2fc9fa9SBarry Smith /* 85c2fc9fa9SBarry Smith DMCoarsen_SNESVI - Computes the regular coarsened DM then computes additional information about its inactive set 86c2fc9fa9SBarry Smith */ 879371c9d4SSatish Balay PetscErrorCode DMCoarsen_SNESVI(DM dm1, MPI_Comm comm, DM *dm2) { 88c2fc9fa9SBarry Smith PetscContainer isnes; 89c2fc9fa9SBarry Smith DM_SNESVI *dmsnesvi1; 90c2fc9fa9SBarry Smith Vec finemarked, coarsemarked; 91c2fc9fa9SBarry Smith IS inactive; 926dbf9973SLawrence Mitchell Mat inject; 93c2fc9fa9SBarry Smith const PetscInt *index; 94c2fc9fa9SBarry Smith PetscInt n, k, cnt = 0, rstart, *coarseindex; 95c2fc9fa9SBarry Smith PetscScalar *marked; 96c2fc9fa9SBarry Smith 97c2fc9fa9SBarry Smith PetscFunctionBegin; 989566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm1, "VI", (PetscObject *)&isnes)); 9928b400f6SJacob Faibussowitsch PetscCheck(isnes, PetscObjectComm((PetscObject)dm1), PETSC_ERR_PLIB, "Composed VI data structure is missing"); 1009566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(isnes, (void **)&dmsnesvi1)); 101c2fc9fa9SBarry Smith 102c2fc9fa9SBarry Smith /* get the original coarsen */ 1039566063dSJacob Faibussowitsch PetscCall((*dmsnesvi1->coarsen)(dm1, comm, dm2)); 104c2fc9fa9SBarry Smith 105c2fc9fa9SBarry Smith /* not sure why this extra reference is needed, but without the dm2 disappears too early */ 10694c98981SBarry Smith /* Updating the KSPCreateVecs() to avoid using DMGetGlobalVector() when matrix is available removes the need for this reference? */ 1079566063dSJacob Faibussowitsch /* PetscCall(PetscObjectReference((PetscObject)*dm2));*/ 108c2fc9fa9SBarry Smith 109c2fc9fa9SBarry Smith /* need to set back global vectors in order to use the original injection */ 1109566063dSJacob Faibussowitsch PetscCall(DMClearGlobalVectors(dm1)); 1111aa26658SKarl Rupp 112c2fc9fa9SBarry Smith dm1->ops->createglobalvector = dmsnesvi1->createglobalvector; 1131aa26658SKarl Rupp 1149566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm1, &finemarked)); 1159566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(*dm2, &coarsemarked)); 116c2fc9fa9SBarry Smith 117c2fc9fa9SBarry Smith /* 118c2fc9fa9SBarry Smith fill finemarked with locations of inactive points 119c2fc9fa9SBarry Smith */ 1209566063dSJacob Faibussowitsch PetscCall(ISGetIndices(dmsnesvi1->inactive, &index)); 1219566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(dmsnesvi1->inactive, &n)); 1229566063dSJacob Faibussowitsch PetscCall(VecSet(finemarked, 0.0)); 12348a46eb9SPierre Jolivet for (k = 0; k < n; k++) PetscCall(VecSetValue(finemarked, index[k], 1.0, INSERT_VALUES)); 1249566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(finemarked)); 1259566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(finemarked)); 126c2fc9fa9SBarry Smith 1279566063dSJacob Faibussowitsch PetscCall(DMCreateInjection(*dm2, dm1, &inject)); 1289566063dSJacob Faibussowitsch PetscCall(MatRestrict(inject, finemarked, coarsemarked)); 1299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&inject)); 130c2fc9fa9SBarry Smith 131c2fc9fa9SBarry Smith /* 132c2fc9fa9SBarry Smith create index set list of coarse inactive points from coarsemarked 133c2fc9fa9SBarry Smith */ 1349566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coarsemarked, &n)); 1359566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(coarsemarked, &rstart, NULL)); 1369566063dSJacob Faibussowitsch PetscCall(VecGetArray(coarsemarked, &marked)); 137c2fc9fa9SBarry Smith for (k = 0; k < n; k++) { 138c2fc9fa9SBarry Smith if (marked[k] != 0.0) cnt++; 139c2fc9fa9SBarry Smith } 1409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cnt, &coarseindex)); 141c2fc9fa9SBarry Smith cnt = 0; 142c2fc9fa9SBarry Smith for (k = 0; k < n; k++) { 143c2fc9fa9SBarry Smith if (marked[k] != 0.0) coarseindex[cnt++] = k + rstart; 144c2fc9fa9SBarry Smith } 1459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coarsemarked, &marked)); 1469566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)coarsemarked), cnt, coarseindex, PETSC_OWN_POINTER, &inactive)); 147c2fc9fa9SBarry Smith 1489566063dSJacob Faibussowitsch PetscCall(DMClearGlobalVectors(dm1)); 1491aa26658SKarl Rupp 150c2fc9fa9SBarry Smith dm1->ops->createglobalvector = DMCreateGlobalVector_SNESVI; 1511aa26658SKarl Rupp 1529566063dSJacob Faibussowitsch PetscCall(DMSetVI(*dm2, inactive)); 153c2fc9fa9SBarry Smith 1549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&finemarked)); 1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coarsemarked)); 1569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&inactive)); 157c2fc9fa9SBarry Smith PetscFunctionReturn(0); 158c2fc9fa9SBarry Smith } 159c2fc9fa9SBarry Smith 1609371c9d4SSatish Balay PetscErrorCode DMDestroy_SNESVI(DM_SNESVI *dmsnesvi) { 161c2fc9fa9SBarry Smith PetscFunctionBegin; 162c2fc9fa9SBarry Smith /* reset the base methods in the DM object that were changed when the DM_SNESVI was reset */ 16325296bd5SBarry Smith dmsnesvi->dm->ops->createinterpolation = dmsnesvi->createinterpolation; 164c2fc9fa9SBarry Smith dmsnesvi->dm->ops->coarsen = dmsnesvi->coarsen; 165c2fc9fa9SBarry Smith dmsnesvi->dm->ops->createglobalvector = dmsnesvi->createglobalvector; 1665a84ad33SLisandro Dalcin dmsnesvi->dm->ops->createinjection = dmsnesvi->createinjection; 1674a7a4c06SLawrence Mitchell dmsnesvi->dm->ops->hascreateinjection = dmsnesvi->hascreateinjection; 168c2fc9fa9SBarry Smith /* need to clear out this vectors because some of them may not have a reference to the DM 169c2fc9fa9SBarry Smith but they are counted as having references to the DM in DMDestroy() */ 1709566063dSJacob Faibussowitsch PetscCall(DMClearGlobalVectors(dmsnesvi->dm)); 171c2fc9fa9SBarry Smith 1729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dmsnesvi->inactive)); 1739566063dSJacob Faibussowitsch PetscCall(PetscFree(dmsnesvi)); 174c2fc9fa9SBarry Smith PetscFunctionReturn(0); 175c2fc9fa9SBarry Smith } 176c2fc9fa9SBarry Smith 177*f6dfbefdSBarry Smith /*@ 178c2fc9fa9SBarry Smith DMSetVI - Marks a DM as associated with a VI problem. This causes the interpolation/restriction operators to 179c2fc9fa9SBarry Smith be restricted to only those variables NOT associated with active constraints. 180c2fc9fa9SBarry Smith 181*f6dfbefdSBarry Smith Logically Collective on dm 182*f6dfbefdSBarry Smith 183*f6dfbefdSBarry Smith Input Parameters: 184*f6dfbefdSBarry Smith + dm - the `DM` object 185*f6dfbefdSBarry Smith - inactive - an `IS` indicating which points are currently not active 186*f6dfbefdSBarry Smith 187*f6dfbefdSBarry Smith .seealso: `SNESVINEWTONRSLS`, `SNESVIGetInactiveSet()` 188*f6dfbefdSBarry Smith @*/ 189*f6dfbefdSBarry Smith PetscErrorCode DMSetVI(DM dm, IS inactive) { 190c2fc9fa9SBarry Smith PetscContainer isnes; 191c2fc9fa9SBarry Smith DM_SNESVI *dmsnesvi; 192c2fc9fa9SBarry Smith 193c2fc9fa9SBarry Smith PetscFunctionBegin; 194c2fc9fa9SBarry Smith if (!dm) PetscFunctionReturn(0); 195c2fc9fa9SBarry Smith 1969566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)inactive)); 197c2fc9fa9SBarry Smith 1989566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "VI", (PetscObject *)&isnes)); 199c2fc9fa9SBarry Smith if (!isnes) { 2009566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)dm), &isnes)); 2019566063dSJacob Faibussowitsch PetscCall(PetscContainerSetUserDestroy(isnes, (PetscErrorCode(*)(void *))DMDestroy_SNESVI)); 2029566063dSJacob Faibussowitsch PetscCall(PetscNew(&dmsnesvi)); 2039566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(isnes, (void *)dmsnesvi)); 2049566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "VI", (PetscObject)isnes)); 2059566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&isnes)); 2061aa26658SKarl Rupp 20725296bd5SBarry Smith dmsnesvi->createinterpolation = dm->ops->createinterpolation; 20825296bd5SBarry Smith dm->ops->createinterpolation = DMCreateInterpolation_SNESVI; 209c2fc9fa9SBarry Smith dmsnesvi->coarsen = dm->ops->coarsen; 210c2fc9fa9SBarry Smith dm->ops->coarsen = DMCoarsen_SNESVI; 211c2fc9fa9SBarry Smith dmsnesvi->createglobalvector = dm->ops->createglobalvector; 212c2fc9fa9SBarry Smith dm->ops->createglobalvector = DMCreateGlobalVector_SNESVI; 2135a84ad33SLisandro Dalcin dmsnesvi->createinjection = dm->ops->createinjection; 2145a84ad33SLisandro Dalcin dm->ops->createinjection = NULL; 2154a7a4c06SLawrence Mitchell dmsnesvi->hascreateinjection = dm->ops->hascreateinjection; 2165a84ad33SLisandro Dalcin dm->ops->hascreateinjection = NULL; 217c2fc9fa9SBarry Smith } else { 2189566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(isnes, (void **)&dmsnesvi)); 2199566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dmsnesvi->inactive)); 220c2fc9fa9SBarry Smith } 2219566063dSJacob Faibussowitsch PetscCall(DMClearGlobalVectors(dm)); 2229566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(inactive, &dmsnesvi->n)); 2231aa26658SKarl Rupp 224c2fc9fa9SBarry Smith dmsnesvi->inactive = inactive; 225c2fc9fa9SBarry Smith dmsnesvi->dm = dm; 226c2fc9fa9SBarry Smith PetscFunctionReturn(0); 227c2fc9fa9SBarry Smith } 228c2fc9fa9SBarry Smith 229c2fc9fa9SBarry Smith /* 230c2fc9fa9SBarry Smith DMDestroyVI - Frees the DM_SNESVI object contained in the DM 23125296bd5SBarry Smith - also resets the function pointers in the DM for createinterpolation() etc to use the original DM 232c2fc9fa9SBarry Smith */ 233*f6dfbefdSBarry Smith PetscErrorCode DMDestroyVI(DM dm) { 234c2fc9fa9SBarry Smith PetscFunctionBegin; 235c2fc9fa9SBarry Smith if (!dm) PetscFunctionReturn(0); 2369566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "VI", (PetscObject)NULL)); 237c2fc9fa9SBarry Smith PetscFunctionReturn(0); 238c2fc9fa9SBarry Smith } 239c2fc9fa9SBarry Smith 2409371c9d4SSatish Balay PetscErrorCode SNESCreateIndexSets_VINEWTONRSLS(SNES snes, Vec X, Vec F, IS *ISact, IS *ISinact) { 241c2fc9fa9SBarry Smith PetscFunctionBegin; 2429566063dSJacob Faibussowitsch PetscCall(SNESVIGetActiveSetIS(snes, X, F, ISact)); 2439566063dSJacob Faibussowitsch PetscCall(ISComplement(*ISact, X->map->rstart, X->map->rend, ISinact)); 244c2fc9fa9SBarry Smith PetscFunctionReturn(0); 245c2fc9fa9SBarry Smith } 246c2fc9fa9SBarry Smith 247c2fc9fa9SBarry Smith /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */ 2489371c9d4SSatish Balay PetscErrorCode SNESCreateSubVectors_VINEWTONRSLS(SNES snes, PetscInt n, Vec *newv) { 249c2fc9fa9SBarry Smith Vec v; 250c2fc9fa9SBarry Smith 251c2fc9fa9SBarry Smith PetscFunctionBegin; 2529566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)snes), &v)); 2539566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v, n, PETSC_DECIDE)); 2549566063dSJacob Faibussowitsch PetscCall(VecSetType(v, VECSTANDARD)); 255c2fc9fa9SBarry Smith *newv = v; 256c2fc9fa9SBarry Smith PetscFunctionReturn(0); 257c2fc9fa9SBarry Smith } 258c2fc9fa9SBarry Smith 259c2fc9fa9SBarry Smith /* Resets the snes PC and KSP when the active set sizes change */ 2609371c9d4SSatish Balay PetscErrorCode SNESVIResetPCandKSP(SNES snes, Mat Amat, Mat Pmat) { 261c2fc9fa9SBarry Smith KSP snesksp; 262c2fc9fa9SBarry Smith 263c2fc9fa9SBarry Smith PetscFunctionBegin; 2649566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &snesksp)); 2659566063dSJacob Faibussowitsch PetscCall(KSPReset(snesksp)); 2669566063dSJacob Faibussowitsch PetscCall(KSPResetFromOptions(snesksp)); 267c2fc9fa9SBarry Smith 268c2fc9fa9SBarry Smith /* 269c2fc9fa9SBarry Smith KSP kspnew; 270c2fc9fa9SBarry Smith PC pcnew; 271ea799195SBarry Smith MatSolverType stype; 272c2fc9fa9SBarry Smith 2739566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)snes),&kspnew)); 274c2fc9fa9SBarry Smith kspnew->pc_side = snesksp->pc_side; 275c2fc9fa9SBarry Smith kspnew->rtol = snesksp->rtol; 276c2fc9fa9SBarry Smith kspnew->abstol = snesksp->abstol; 277c2fc9fa9SBarry Smith kspnew->max_it = snesksp->max_it; 2789566063dSJacob Faibussowitsch PetscCall(KSPSetType(kspnew,((PetscObject)snesksp)->type_name)); 2799566063dSJacob Faibussowitsch PetscCall(KSPGetPC(kspnew,&pcnew)); 2809566063dSJacob Faibussowitsch PetscCall(PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name)); 2819566063dSJacob Faibussowitsch PetscCall(PCSetOperators(kspnew->pc,Amat,Pmat)); 2829566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(snesksp->pc,&stype)); 2839566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(kspnew->pc,stype)); 2849566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&snesksp)); 285c2fc9fa9SBarry Smith snes->ksp = kspnew; 2869566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)snes,(PetscObject)kspnew)); 2879566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(kspnew));*/ 288c2fc9fa9SBarry Smith PetscFunctionReturn(0); 289c2fc9fa9SBarry Smith } 290c2fc9fa9SBarry Smith 291c2fc9fa9SBarry Smith /* Variational Inequality solver using reduce space method. No semismooth algorithm is 292c2fc9fa9SBarry Smith implemented in this algorithm. It basically identifies the active constraints and does 293c2fc9fa9SBarry Smith a linear solve on the other variables (those not associated with the active constraints). */ 2949371c9d4SSatish Balay PetscErrorCode SNESSolve_VINEWTONRSLS(SNES snes) { 295f450aa47SBarry Smith SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data; 296c2fc9fa9SBarry Smith PetscInt maxits, i, lits; 297422a814eSBarry Smith SNESLineSearchReason lssucceed; 298c2fc9fa9SBarry Smith PetscReal fnorm, gnorm, xnorm = 0, ynorm; 2999bd66eb0SPeter Brune Vec Y, X, F; 300c2fc9fa9SBarry Smith KSPConvergedReason kspreason; 30192e89061SBarry Smith KSP ksp; 30292e89061SBarry Smith PC pc; 303c2fc9fa9SBarry Smith 304c2fc9fa9SBarry Smith PetscFunctionBegin; 30592e89061SBarry Smith /* Multigrid must use Galerkin for coarse grids with active set/reduced space methods; cannot rediscretize on coarser grids*/ 3069566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 3079566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 3089566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_BOTH)); 30992e89061SBarry Smith 310c2fc9fa9SBarry Smith snes->numFailures = 0; 311c2fc9fa9SBarry Smith snes->numLinearSolveFailures = 0; 312c2fc9fa9SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 313c2fc9fa9SBarry Smith 314c2fc9fa9SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 315c2fc9fa9SBarry Smith X = snes->vec_sol; /* solution vector */ 316c2fc9fa9SBarry Smith F = snes->vec_func; /* residual vector */ 317c2fc9fa9SBarry Smith Y = snes->work[0]; /* work vectors */ 3189bd66eb0SPeter Brune 3199566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetVIFunctions(snes->linesearch, SNESVIProjectOntoBounds, SNESVIComputeInactiveSetFnorm)); 3209566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetVecs(snes->linesearch, X, NULL, NULL, NULL, NULL)); 3219566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetUp(snes->linesearch)); 322c2fc9fa9SBarry Smith 3239566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 324c2fc9fa9SBarry Smith snes->iter = 0; 325c2fc9fa9SBarry Smith snes->norm = 0.0; 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 327c2fc9fa9SBarry Smith 3289566063dSJacob Faibussowitsch PetscCall(SNESVIProjectOntoBounds(snes, X)); 3299566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 3309566063dSJacob Faibussowitsch PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm)); 3319566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- ||x|| */ 332422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 3339566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 334c2fc9fa9SBarry Smith snes->norm = fnorm; 3359566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3369566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 3379566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 338c2fc9fa9SBarry Smith 339c2fc9fa9SBarry Smith /* test convergence */ 340dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 341c2fc9fa9SBarry Smith if (snes->reason) PetscFunctionReturn(0); 342c2fc9fa9SBarry Smith 343c2fc9fa9SBarry Smith for (i = 0; i < maxits; i++) { 344f009fc93SPatrick Farrell IS IS_act; /* _act -> active set _inact -> inactive set */ 345c2fc9fa9SBarry Smith IS IS_redact; /* redundant active set */ 346c2fc9fa9SBarry Smith VecScatter scat_act, scat_inact; 347c2fc9fa9SBarry Smith PetscInt nis_act, nis_inact; 348c2fc9fa9SBarry Smith Vec Y_act, Y_inact, F_inact; 349c2fc9fa9SBarry Smith Mat jac_inact_inact, prejac_inact_inact; 350c2fc9fa9SBarry Smith PetscBool isequal; 351c2fc9fa9SBarry Smith 352c2fc9fa9SBarry Smith /* Call general purpose update function */ 353dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 3549566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 35507b62357SFande Kong SNESCheckJacobianDomainerror(snes); 356c2fc9fa9SBarry Smith 357c2fc9fa9SBarry Smith /* Create active and inactive index sets */ 358c2fc9fa9SBarry Smith 359c2fc9fa9SBarry Smith /*original 3609566063dSJacob Faibussowitsch PetscCall(SNESVICreateIndexSets_RS(snes,X,F,&IS_act,&vi->IS_inact)); 361c2fc9fa9SBarry Smith */ 3629566063dSJacob Faibussowitsch PetscCall(SNESVIGetActiveSetIS(snes, X, F, &IS_act)); 363c2fc9fa9SBarry Smith 364c2fc9fa9SBarry Smith if (vi->checkredundancy) { 3659566063dSJacob Faibussowitsch PetscCall((*vi->checkredundancy)(snes, IS_act, &IS_redact, vi->ctxP)); 366c2fc9fa9SBarry Smith if (IS_redact) { 3679566063dSJacob Faibussowitsch PetscCall(ISSort(IS_redact)); 3689566063dSJacob Faibussowitsch PetscCall(ISComplement(IS_redact, X->map->rstart, X->map->rend, &vi->IS_inact)); 3699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&IS_redact)); 3701aa26658SKarl Rupp } else { 3719566063dSJacob Faibussowitsch PetscCall(ISComplement(IS_act, X->map->rstart, X->map->rend, &vi->IS_inact)); 372c2fc9fa9SBarry Smith } 373c2fc9fa9SBarry Smith } else { 3749566063dSJacob Faibussowitsch PetscCall(ISComplement(IS_act, X->map->rstart, X->map->rend, &vi->IS_inact)); 375c2fc9fa9SBarry Smith } 376c2fc9fa9SBarry Smith 377c2fc9fa9SBarry Smith /* Create inactive set submatrix */ 3789566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(snes->jacobian, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &jac_inact_inact)); 379c2fc9fa9SBarry Smith 38066bfb381SJed Brown if (0) { /* Dead code (temporary developer hack) */ 38166bfb381SJed Brown IS keptrows; 3829566063dSJacob Faibussowitsch PetscCall(MatFindNonzeroRows(jac_inact_inact, &keptrows)); 38366bfb381SJed Brown if (keptrows) { 384c2fc9fa9SBarry Smith PetscInt cnt, *nrows, k; 385c2fc9fa9SBarry Smith const PetscInt *krows, *inact; 386367daffbSBarry Smith PetscInt rstart; 387c2fc9fa9SBarry Smith 3889566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(jac_inact_inact, &rstart, NULL)); 3899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac_inact_inact)); 3909566063dSJacob Faibussowitsch PetscCall(ISDestroy(&IS_act)); 391c2fc9fa9SBarry Smith 3929566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(keptrows, &cnt)); 3939566063dSJacob Faibussowitsch PetscCall(ISGetIndices(keptrows, &krows)); 3949566063dSJacob Faibussowitsch PetscCall(ISGetIndices(vi->IS_inact, &inact)); 3959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cnt, &nrows)); 3961aa26658SKarl Rupp for (k = 0; k < cnt; k++) nrows[k] = inact[krows[k] - rstart]; 3979566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(keptrows, &krows)); 3989566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(vi->IS_inact, &inact)); 3999566063dSJacob Faibussowitsch PetscCall(ISDestroy(&keptrows)); 4009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vi->IS_inact)); 401c2fc9fa9SBarry Smith 4029566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), cnt, nrows, PETSC_OWN_POINTER, &vi->IS_inact)); 4039566063dSJacob Faibussowitsch PetscCall(ISComplement(vi->IS_inact, F->map->rstart, F->map->rend, &IS_act)); 4049566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(snes->jacobian, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &jac_inact_inact)); 405c2fc9fa9SBarry Smith } 40666bfb381SJed Brown } 4079566063dSJacob Faibussowitsch PetscCall(DMSetVI(snes->dm, vi->IS_inact)); 408c2fc9fa9SBarry Smith /* remove later */ 409c2fc9fa9SBarry Smith 410c2fc9fa9SBarry Smith /* 4119566063dSJacob Faibussowitsch PetscCall(VecView(vi->xu,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xu))->comm))); 4129566063dSJacob Faibussowitsch PetscCall(VecView(vi->xl,PETSC_VIEWER_BINARY_(((PetscObject)(vi->xl))->comm))); 4139566063dSJacob Faibussowitsch PetscCall(VecView(X,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)X)))); 4149566063dSJacob Faibussowitsch PetscCall(VecView(F,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)F)))); 4159566063dSJacob Faibussowitsch PetscCall(ISView(vi->IS_inact,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)vi->IS_inact)))); 416c2fc9fa9SBarry Smith */ 417c2fc9fa9SBarry Smith 418c2fc9fa9SBarry Smith /* Get sizes of active and inactive sets */ 4199566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(IS_act, &nis_act)); 4209566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(vi->IS_inact, &nis_inact)); 421c2fc9fa9SBarry Smith 422c2fc9fa9SBarry Smith /* Create active and inactive set vectors */ 4239566063dSJacob Faibussowitsch PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_inact, &F_inact)); 4249566063dSJacob Faibussowitsch PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_act, &Y_act)); 4259566063dSJacob Faibussowitsch PetscCall(SNESCreateSubVectors_VINEWTONRSLS(snes, nis_inact, &Y_inact)); 426c2fc9fa9SBarry Smith 427c2fc9fa9SBarry Smith /* Create scatter contexts */ 4289566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(Y, IS_act, Y_act, NULL, &scat_act)); 4299566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(Y, vi->IS_inact, Y_inact, NULL, &scat_inact)); 430c2fc9fa9SBarry Smith 431c2fc9fa9SBarry Smith /* Do a vec scatter to active and inactive set vectors */ 4329566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat_inact, F, F_inact, INSERT_VALUES, SCATTER_FORWARD)); 4339566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat_inact, F, F_inact, INSERT_VALUES, SCATTER_FORWARD)); 434c2fc9fa9SBarry Smith 4359566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat_act, Y, Y_act, INSERT_VALUES, SCATTER_FORWARD)); 4369566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat_act, Y, Y_act, INSERT_VALUES, SCATTER_FORWARD)); 437c2fc9fa9SBarry Smith 4389566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat_inact, Y, Y_inact, INSERT_VALUES, SCATTER_FORWARD)); 4399566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat_inact, Y, Y_inact, INSERT_VALUES, SCATTER_FORWARD)); 440c2fc9fa9SBarry Smith 441c2fc9fa9SBarry Smith /* Active set direction = 0 */ 4429566063dSJacob Faibussowitsch PetscCall(VecSet(Y_act, 0)); 443c2fc9fa9SBarry Smith if (snes->jacobian != snes->jacobian_pre) { 4449566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(snes->jacobian_pre, vi->IS_inact, vi->IS_inact, MAT_INITIAL_MATRIX, &prejac_inact_inact)); 445c2fc9fa9SBarry Smith } else prejac_inact_inact = jac_inact_inact; 446c2fc9fa9SBarry Smith 4479566063dSJacob Faibussowitsch PetscCall(ISEqual(vi->IS_inact_prev, vi->IS_inact, &isequal)); 448c2fc9fa9SBarry Smith if (!isequal) { 4499566063dSJacob Faibussowitsch PetscCall(SNESVIResetPCandKSP(snes, jac_inact_inact, prejac_inact_inact)); 4509566063dSJacob Faibussowitsch PetscCall(PCFieldSplitRestrictIS(pc, vi->IS_inact)); 451c2fc9fa9SBarry Smith } 452c2fc9fa9SBarry Smith 4539566063dSJacob Faibussowitsch /* PetscCall(ISView(vi->IS_inact,0)); */ 4549566063dSJacob Faibussowitsch /* PetscCall(ISView(IS_act,0));*/ 455c2fc9fa9SBarry Smith /* ierr = MatView(snes->jacobian_pre,0); */ 456c2fc9fa9SBarry Smith 4579566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, jac_inact_inact, prejac_inact_inact)); 4589566063dSJacob Faibussowitsch PetscCall(KSPSetUp(snes->ksp)); 459c2fc9fa9SBarry Smith { 460c2fc9fa9SBarry Smith PC pc; 461c2fc9fa9SBarry Smith PetscBool flg; 4629566063dSJacob Faibussowitsch PetscCall(KSPGetPC(snes->ksp, &pc)); 4639566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &flg)); 464c2fc9fa9SBarry Smith if (flg) { 465c2fc9fa9SBarry Smith KSP *subksps; 4669566063dSJacob Faibussowitsch PetscCall(PCFieldSplitGetSubKSP(pc, NULL, &subksps)); 4679566063dSJacob Faibussowitsch PetscCall(KSPGetPC(subksps[0], &pc)); 4689566063dSJacob Faibussowitsch PetscCall(PetscFree(subksps)); 4699566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &flg)); 470c2fc9fa9SBarry Smith if (flg) { 471c2fc9fa9SBarry Smith PetscInt n, N = 101 * 101, j, cnts[3] = {0, 0, 0}; 472c2fc9fa9SBarry Smith const PetscInt *ii; 473c2fc9fa9SBarry Smith 4749566063dSJacob Faibussowitsch PetscCall(ISGetSize(vi->IS_inact, &n)); 4759566063dSJacob Faibussowitsch PetscCall(ISGetIndices(vi->IS_inact, &ii)); 476c2fc9fa9SBarry Smith for (j = 0; j < n; j++) { 477c2fc9fa9SBarry Smith if (ii[j] < N) cnts[0]++; 478c2fc9fa9SBarry Smith else if (ii[j] < 2 * N) cnts[1]++; 479c2fc9fa9SBarry Smith else if (ii[j] < 3 * N) cnts[2]++; 480c2fc9fa9SBarry Smith } 4819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(vi->IS_inact, &ii)); 482c2fc9fa9SBarry Smith 4839566063dSJacob Faibussowitsch PetscCall(PCBJacobiSetTotalBlocks(pc, 3, cnts)); 484c2fc9fa9SBarry Smith } 485c2fc9fa9SBarry Smith } 486c2fc9fa9SBarry Smith } 487c2fc9fa9SBarry Smith 4889566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F_inact, Y_inact)); 4899566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat_act, Y_act, Y, INSERT_VALUES, SCATTER_REVERSE)); 4909566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat_act, Y_act, Y, INSERT_VALUES, SCATTER_REVERSE)); 4919566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat_inact, Y_inact, Y, INSERT_VALUES, SCATTER_REVERSE)); 4929566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat_inact, Y_inact, Y, INSERT_VALUES, SCATTER_REVERSE)); 493c2fc9fa9SBarry Smith 4949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F_inact)); 4959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y_act)); 4969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y_inact)); 4979566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scat_act)); 4989566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scat_inact)); 4999566063dSJacob Faibussowitsch PetscCall(ISDestroy(&IS_act)); 500c2fc9fa9SBarry Smith if (!isequal) { 5019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vi->IS_inact_prev)); 5029566063dSJacob Faibussowitsch PetscCall(ISDuplicate(vi->IS_inact, &vi->IS_inact_prev)); 503c2fc9fa9SBarry Smith } 5049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vi->IS_inact)); 5059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac_inact_inact)); 50648a46eb9SPierre Jolivet if (snes->jacobian != snes->jacobian_pre) PetscCall(MatDestroy(&prejac_inact_inact)); 507fa6eefd1SCian Wilson 5089566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason)); 509fa6eefd1SCian Wilson if (kspreason < 0) { 510fa6eefd1SCian Wilson if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 51163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", number linear solve failures %" PetscInt_FMT " greater than current SNES allowed, stopping solve\n", snes->iter, snes->numLinearSolveFailures)); 512fa6eefd1SCian Wilson snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 513fa6eefd1SCian Wilson break; 514fa6eefd1SCian Wilson } 515fa6eefd1SCian Wilson } 516fa6eefd1SCian Wilson 5179566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 518c2fc9fa9SBarry Smith snes->linear_its += lits; 51963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 520c2fc9fa9SBarry Smith /* 5216b2b7091SBarry Smith if (snes->ops->precheck) { 522c2fc9fa9SBarry Smith PetscBool changed_y = PETSC_FALSE; 523dbbe0bcdSBarry Smith PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y); 524c2fc9fa9SBarry Smith } 525c2fc9fa9SBarry Smith 5261baa6e33SBarry Smith if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W)); 527c2fc9fa9SBarry Smith */ 528c2fc9fa9SBarry Smith /* Compute a (scaled) negative update in the line search routine: 529c2fc9fa9SBarry Smith Y <- X - lambda*Y 530c2fc9fa9SBarry Smith and evaluate G = function(Y) (depends on the line search). 531c2fc9fa9SBarry Smith */ 5329566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, snes->vec_sol_update)); 5339371c9d4SSatish Balay ynorm = 1; 5349371c9d4SSatish Balay gnorm = fnorm; 5359566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &gnorm, Y)); 5369566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 5379566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm)); 5389566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lssucceed)); 539c2fc9fa9SBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 540c2fc9fa9SBarry Smith if (snes->domainerror) { 541c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 5429566063dSJacob Faibussowitsch PetscCall(DMDestroyVI(snes->dm)); 543c2fc9fa9SBarry Smith PetscFunctionReturn(0); 544c2fc9fa9SBarry Smith } 545422a814eSBarry Smith if (lssucceed) { 546c2fc9fa9SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 547c2fc9fa9SBarry Smith PetscBool ismin; 548c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 5499566063dSJacob Faibussowitsch PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, F, X, gnorm, &ismin)); 550c2fc9fa9SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 551c2fc9fa9SBarry Smith break; 552c2fc9fa9SBarry Smith } 553c2fc9fa9SBarry Smith } 5549566063dSJacob Faibussowitsch PetscCall(DMDestroyVI(snes->dm)); 555c2fc9fa9SBarry Smith /* Update function and solution vectors */ 556c2fc9fa9SBarry Smith fnorm = gnorm; 557c2fc9fa9SBarry Smith /* Monitor convergence */ 5589566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 559c2fc9fa9SBarry Smith snes->iter = i + 1; 560c2fc9fa9SBarry Smith snes->norm = fnorm; 561c1e67a49SFande Kong snes->xnorm = xnorm; 562c1e67a49SFande Kong snes->ynorm = ynorm; 5639566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5649566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 5659566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 566c2fc9fa9SBarry Smith /* Test for convergence, xnorm = || X || */ 5679566063dSJacob Faibussowitsch if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 568dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 569c2fc9fa9SBarry Smith if (snes->reason) break; 570c2fc9fa9SBarry Smith } 57191a42fcfSBarry Smith /* 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 */ 5729566063dSJacob Faibussowitsch PetscCall(DMDestroyVI(snes->dm)); 573c2fc9fa9SBarry Smith if (i == maxits) { 57463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 575c2fc9fa9SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 576c2fc9fa9SBarry Smith } 577c2fc9fa9SBarry Smith PetscFunctionReturn(0); 578c2fc9fa9SBarry Smith } 579c2fc9fa9SBarry Smith 580*f6dfbefdSBarry Smith /*@C 581*f6dfbefdSBarry Smith SNESVISetRedundancyCheck - Provide a function to check for any redundancy in the VI active set 582*f6dfbefdSBarry Smith 583*f6dfbefdSBarry Smith Logically Collective on snes 584*f6dfbefdSBarry Smith 585*f6dfbefdSBarry Smith Input Parameters: 586*f6dfbefdSBarry Smith + snes - the `SNESVINEWTONRSLS` context 587*f6dfbefdSBarry Smith . func - the function to check of redundancies 588*f6dfbefdSBarry Smith - ctx - optional context used by the function 589*f6dfbefdSBarry Smith 590*f6dfbefdSBarry Smith Level: advanced 591*f6dfbefdSBarry Smith 592*f6dfbefdSBarry Smith Note: 593*f6dfbefdSBarry Smith Sometimes the inactive set will result in a non-singular sub-Jacobian problem that needs to be solved, this allows the user, 594*f6dfbefdSBarry Smith when they know more about their specific problem to provide a function that removes the redundancy that results in the singular linear system 595*f6dfbefdSBarry Smith 596*f6dfbefdSBarry Smith .seealso: `SNESVINEWTONRSLS`, `SNESVIGetInactiveSet()`, `DMSetVI()` 597*f6dfbefdSBarry Smith @*/ 5989371c9d4SSatish Balay PetscErrorCode SNESVISetRedundancyCheck(SNES snes, PetscErrorCode (*func)(SNES, IS, IS *, void *), void *ctx) { 599f450aa47SBarry Smith SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data; 600c2fc9fa9SBarry Smith 601c2fc9fa9SBarry Smith PetscFunctionBegin; 602c2fc9fa9SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 603c2fc9fa9SBarry Smith vi->checkredundancy = func; 604c2fc9fa9SBarry Smith vi->ctxP = ctx; 605c2fc9fa9SBarry Smith PetscFunctionReturn(0); 606c2fc9fa9SBarry Smith } 607c2fc9fa9SBarry Smith 608c2fc9fa9SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 609c2fc9fa9SBarry Smith #include <engine.h> 610c2fc9fa9SBarry Smith #include <mex.h> 6119371c9d4SSatish Balay typedef struct { 6129371c9d4SSatish Balay char *funcname; 6139371c9d4SSatish Balay mxArray *ctx; 6149371c9d4SSatish Balay } SNESMatlabContext; 615c2fc9fa9SBarry Smith 6169371c9d4SSatish Balay PetscErrorCode SNESVIRedundancyCheck_Matlab(SNES snes, IS is_act, IS *is_redact, void *ctx) { 617c2fc9fa9SBarry Smith SNESMatlabContext *sctx = (SNESMatlabContext *)ctx; 618c2fc9fa9SBarry Smith int nlhs = 1, nrhs = 5; 619c2fc9fa9SBarry Smith mxArray *plhs[1], *prhs[5]; 620c2fc9fa9SBarry Smith long long int l1 = 0, l2 = 0, ls = 0; 6210298fd71SBarry Smith PetscInt *indices = NULL; 622c2fc9fa9SBarry Smith 623c2fc9fa9SBarry Smith PetscFunctionBegin; 624c2fc9fa9SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 625c2fc9fa9SBarry Smith PetscValidHeaderSpecific(is_act, IS_CLASSID, 2); 626c2fc9fa9SBarry Smith PetscValidPointer(is_redact, 3); 627c2fc9fa9SBarry Smith PetscCheckSameComm(snes, 1, is_act, 2); 628c2fc9fa9SBarry Smith 629c2fc9fa9SBarry Smith /* Create IS for reduced active set of size 0, its size and indices will 630c2fc9fa9SBarry Smith bet set by the Matlab function */ 6319566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), 0, indices, PETSC_OWN_POINTER, is_redact)); 632c2fc9fa9SBarry Smith /* call Matlab function in ctx */ 6339566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&ls, &snes, 1)); 6349566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&l1, &is_act, 1)); 6359566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&l2, is_redact, 1)); 636c2fc9fa9SBarry Smith prhs[0] = mxCreateDoubleScalar((double)ls); 637c2fc9fa9SBarry Smith prhs[1] = mxCreateDoubleScalar((double)l1); 638c2fc9fa9SBarry Smith prhs[2] = mxCreateDoubleScalar((double)l2); 639c2fc9fa9SBarry Smith prhs[3] = mxCreateString(sctx->funcname); 640c2fc9fa9SBarry Smith prhs[4] = sctx->ctx; 6419566063dSJacob Faibussowitsch PetscCall(mexCallMATLAB(nlhs, plhs, nrhs, prhs, "PetscSNESVIRedundancyCheckInternal")); 6429566063dSJacob Faibussowitsch PetscCall(mxGetScalar(plhs[0])); 643c2fc9fa9SBarry Smith mxDestroyArray(prhs[0]); 644c2fc9fa9SBarry Smith mxDestroyArray(prhs[1]); 645c2fc9fa9SBarry Smith mxDestroyArray(prhs[2]); 646c2fc9fa9SBarry Smith mxDestroyArray(prhs[3]); 647c2fc9fa9SBarry Smith mxDestroyArray(plhs[0]); 648c2fc9fa9SBarry Smith PetscFunctionReturn(0); 649c2fc9fa9SBarry Smith } 650c2fc9fa9SBarry Smith 6519371c9d4SSatish Balay PetscErrorCode SNESVISetRedundancyCheckMatlab(SNES snes, const char *func, mxArray *ctx) { 652c2fc9fa9SBarry Smith SNESMatlabContext *sctx; 653c2fc9fa9SBarry Smith 654c2fc9fa9SBarry Smith PetscFunctionBegin; 655c2fc9fa9SBarry Smith /* currently sctx is memory bleed */ 6569566063dSJacob Faibussowitsch PetscCall(PetscNew(&sctx)); 6579566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(func, &sctx->funcname)); 658c2fc9fa9SBarry Smith sctx->ctx = mxDuplicateArray(ctx); 6599566063dSJacob Faibussowitsch PetscCall(SNESVISetRedundancyCheck(snes, SNESVIRedundancyCheck_Matlab, sctx)); 660c2fc9fa9SBarry Smith PetscFunctionReturn(0); 661c2fc9fa9SBarry Smith } 662c2fc9fa9SBarry Smith 663c2fc9fa9SBarry Smith #endif 664c2fc9fa9SBarry Smith 665c2fc9fa9SBarry Smith /* 666f450aa47SBarry Smith SNESSetUp_VINEWTONRSLS - Sets up the internal data structures for the later use 667c2fc9fa9SBarry Smith of the SNESVI nonlinear solver. 668c2fc9fa9SBarry Smith 669c2fc9fa9SBarry Smith Input Parameter: 670c2fc9fa9SBarry Smith . snes - the SNES context 671c2fc9fa9SBarry Smith 672c2fc9fa9SBarry Smith Application Interface Routine: SNESSetUp() 673c2fc9fa9SBarry Smith 674*f6dfbefdSBarry Smith Note: 675c2fc9fa9SBarry Smith For basic use of the SNES solvers, the user need not explicitly call 676c2fc9fa9SBarry Smith SNESSetUp(), since these actions will automatically occur during 677c2fc9fa9SBarry Smith the call to SNESSolve(). 678c2fc9fa9SBarry Smith */ 6799371c9d4SSatish Balay PetscErrorCode SNESSetUp_VINEWTONRSLS(SNES snes) { 680f450aa47SBarry Smith SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data; 681c2fc9fa9SBarry Smith PetscInt *indices; 682c2fc9fa9SBarry Smith PetscInt i, n, rstart, rend; 683f1c6b773SPeter Brune SNESLineSearch linesearch; 684c2fc9fa9SBarry Smith 685c2fc9fa9SBarry Smith PetscFunctionBegin; 6869566063dSJacob Faibussowitsch PetscCall(SNESSetUp_VI(snes)); 687c2fc9fa9SBarry Smith 688c2fc9fa9SBarry Smith /* Set up previous active index set for the first snes solve 689c2fc9fa9SBarry Smith vi->IS_inact_prev = 0,1,2,....N */ 690c2fc9fa9SBarry Smith 6919566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(snes->vec_sol, &rstart, &rend)); 6929566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 6939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &indices)); 694c2fc9fa9SBarry Smith for (i = 0; i < n; i++) indices[i] = rstart + i; 6959566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), n, indices, PETSC_OWN_POINTER, &vi->IS_inact_prev)); 6969bd66eb0SPeter Brune 6979bd66eb0SPeter Brune /* set the line search functions */ 6989bd66eb0SPeter Brune if (!snes->linesearch) { 6999566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 7009566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 7019bd66eb0SPeter Brune } 702c2fc9fa9SBarry Smith PetscFunctionReturn(0); 703c2fc9fa9SBarry Smith } 704*f6dfbefdSBarry Smith 7059371c9d4SSatish Balay PetscErrorCode SNESReset_VINEWTONRSLS(SNES snes) { 706f450aa47SBarry Smith SNES_VINEWTONRSLS *vi = (SNES_VINEWTONRSLS *)snes->data; 707c2fc9fa9SBarry Smith 708c2fc9fa9SBarry Smith PetscFunctionBegin; 7099566063dSJacob Faibussowitsch PetscCall(SNESReset_VI(snes)); 7109566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vi->IS_inact_prev)); 711c2fc9fa9SBarry Smith PetscFunctionReturn(0); 712c2fc9fa9SBarry Smith } 713c2fc9fa9SBarry Smith 714c2fc9fa9SBarry Smith /*MC 715f450aa47SBarry Smith SNESVINEWTONRSLS - Reduced space active set solvers for variational inequalities based on Newton's method 716c2fc9fa9SBarry Smith 717*f6dfbefdSBarry Smith Options Database Keys: 718b621fa8fSRichard Tran Mills + -snes_type <vinewtonssls,vinewtonrsls> - a semi-smooth solver, a reduced space active set method 71961589011SJed Brown - -snes_vi_monitor - prints the number of active constraints at each iteration. 720c2fc9fa9SBarry Smith 721c2fc9fa9SBarry Smith Level: beginner 722c2fc9fa9SBarry Smith 723b80f3ac1SShri Abhyankar References: 724606c0280SSatish Balay . * - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale 725b80f3ac1SShri Abhyankar Applications, Optimization Methods and Software, 21 (2006). 726b80f3ac1SShri Abhyankar 727*f6dfbefdSBarry Smith Note: 728*f6dfbefdSBarry Smith At each set of this methods the algorithm produces an inactive set of variables that are constrained to their current values 729*f6dfbefdSBarry Smith (because changing these values would result in those variables no longer satisfying the inequality constraints) 730*f6dfbefdSBarry Smith and produces a step direction by solving the linear system arising from the Jacobian with the inactive variables removed. In other 731*f6dfbefdSBarry Smith words on a reduced space of the solution space. Based on the Newton update it then adjusts the inactive sep for the next iteration. 732c2fc9fa9SBarry Smith 733*f6dfbefdSBarry Smith .seealso: `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONSSLS`, `SNESNEWTONTRDC`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESVIGetInactiveSet()`, `DMSetVI()`, `SNESVISetRedundancyCheck()` 734c2fc9fa9SBarry Smith M*/ 7359371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONRSLS(SNES snes) { 736f450aa47SBarry Smith SNES_VINEWTONRSLS *vi; 737d8d34be6SBarry Smith SNESLineSearch linesearch; 738c2fc9fa9SBarry Smith 739c2fc9fa9SBarry Smith PetscFunctionBegin; 740f450aa47SBarry Smith snes->ops->reset = SNESReset_VINEWTONRSLS; 741f450aa47SBarry Smith snes->ops->setup = SNESSetUp_VINEWTONRSLS; 742f450aa47SBarry Smith snes->ops->solve = SNESSolve_VINEWTONRSLS; 743c2fc9fa9SBarry Smith snes->ops->destroy = SNESDestroy_VI; 744c2fc9fa9SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_VI; 7450298fd71SBarry Smith snes->ops->view = NULL; 7468d359177SBarry Smith snes->ops->converged = SNESConvergedDefault_VI; 747c2fc9fa9SBarry Smith 748c2fc9fa9SBarry Smith snes->usesksp = PETSC_TRUE; 749efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 750c2fc9fa9SBarry Smith 7519566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 75248a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 7539566063dSJacob Faibussowitsch PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0)); 754d8d34be6SBarry Smith 7554fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 7564fc747eaSLawrence Mitchell 7579566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes, &vi)); 758c2fc9fa9SBarry Smith snes->data = (void *)vi; 7590298fd71SBarry Smith vi->checkredundancy = NULL; 760c2fc9fa9SBarry Smith 7619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI)); 7629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI)); 763c2fc9fa9SBarry Smith PetscFunctionReturn(0); 764c2fc9fa9SBarry Smith } 765