1907a3e9cSStefano Zampini #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2907a3e9cSStefano Zampini #include <petsc/private/dmlabelimpl.h> 3907a3e9cSStefano Zampini #include <petsc/private/sectionimpl.h> 4907a3e9cSStefano Zampini #include <petsc/private/sfimpl.h> 5907a3e9cSStefano Zampini 6e4d5475eSStefano Zampini static PetscErrorCode DMTransferMaterialParameters(DM dm, PetscSF sf, DM odm) 7e4d5475eSStefano Zampini { 8e4d5475eSStefano Zampini Vec A; 9e4d5475eSStefano Zampini 10e4d5475eSStefano Zampini PetscFunctionBegin; 11e4d5475eSStefano Zampini /* TODO handle regions? */ 12e4d5475eSStefano Zampini PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &A)); 13e4d5475eSStefano Zampini if (A) { 14e4d5475eSStefano Zampini DM dmAux, ocdm, odmAux; 15e4d5475eSStefano Zampini Vec oA, oAt; 16e4d5475eSStefano Zampini PetscSection sec, osec; 17e4d5475eSStefano Zampini 18e4d5475eSStefano Zampini PetscCall(VecGetDM(A, &dmAux)); 19e4d5475eSStefano Zampini PetscCall(DMClone(odm, &odmAux)); 20e4d5475eSStefano Zampini PetscCall(DMGetCoordinateDM(odm, &ocdm)); 21e4d5475eSStefano Zampini PetscCall(DMSetCoordinateDM(odmAux, ocdm)); 22e4d5475eSStefano Zampini PetscCall(DMCopyDisc(dmAux, odmAux)); 23e4d5475eSStefano Zampini 24e4d5475eSStefano Zampini PetscCall(DMGetLocalSection(dmAux, &sec)); 25e4d5475eSStefano Zampini PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)sec), &osec)); 26e4d5475eSStefano Zampini PetscCall(VecCreate(PetscObjectComm((PetscObject)A), &oAt)); 27e4d5475eSStefano Zampini PetscCall(DMPlexDistributeField(dmAux, sf, sec, A, osec, oAt)); 28e4d5475eSStefano Zampini PetscCall(DMSetLocalSection(odmAux, osec)); 29e4d5475eSStefano Zampini PetscCall(PetscSectionDestroy(&osec)); 30e4d5475eSStefano Zampini PetscCall(DMCreateLocalVector(odmAux, &oA)); 31e4d5475eSStefano Zampini PetscCall(DMDestroy(&odmAux)); 32e4d5475eSStefano Zampini PetscCall(VecCopy(oAt, oA)); 33e4d5475eSStefano Zampini PetscCall(VecDestroy(&oAt)); 34e4d5475eSStefano Zampini PetscCall(DMSetAuxiliaryVec(odm, NULL, 0, 0, oA)); 35e4d5475eSStefano Zampini PetscCall(VecDestroy(&oA)); 36e4d5475eSStefano Zampini } 37e4d5475eSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 38e4d5475eSStefano Zampini } 39e4d5475eSStefano Zampini 40907a3e9cSStefano Zampini PetscErrorCode DMCreateDomainDecomposition_Plex(DM dm, PetscInt *nsub, char ***names, IS **innerises, IS **outerises, DM **dms) 41907a3e9cSStefano Zampini { 42907a3e9cSStefano Zampini DM odm; 43907a3e9cSStefano Zampini PetscSF migrationSF = NULL, sectionSF; 44907a3e9cSStefano Zampini PetscSection sec, tsec, ogsec, olsec; 45907a3e9cSStefano Zampini PetscInt n, mh, ddovl = 0, pStart, pEnd, ni, no, nl; 46907a3e9cSStefano Zampini PetscDS ds; 47907a3e9cSStefano Zampini DMLabel label; 48907a3e9cSStefano Zampini const char *oname = "__internal_plex_dd_ovl_"; 49907a3e9cSStefano Zampini IS gi_is, li_is, go_is, gl_is, ll_is; 50907a3e9cSStefano Zampini IS gis, lis; 51907a3e9cSStefano Zampini PetscInt rst, ren, c, *gidxs, *lidxs, *tidxs; 52907a3e9cSStefano Zampini Vec gvec; 53907a3e9cSStefano Zampini 54907a3e9cSStefano Zampini PetscFunctionBegin; 55907a3e9cSStefano Zampini n = 1; 56907a3e9cSStefano Zampini if (nsub) *nsub = n; 57907a3e9cSStefano Zampini if (names) PetscCall(PetscCalloc1(n, names)); 58907a3e9cSStefano Zampini if (innerises) PetscCall(PetscCalloc1(n, innerises)); 59907a3e9cSStefano Zampini if (outerises) PetscCall(PetscCalloc1(n, outerises)); 60907a3e9cSStefano Zampini if (dms) PetscCall(PetscCalloc1(n, dms)); 61907a3e9cSStefano Zampini 62907a3e9cSStefano Zampini PetscObjectOptionsBegin((PetscObject)dm); 63907a3e9cSStefano Zampini PetscCall(PetscOptionsBoundedInt("-dm_plex_dd_overlap", "The size of the overlap halo for the subdomains", "DMCreateDomainDecomposition", ddovl, &ddovl, NULL, 0)); 64907a3e9cSStefano Zampini PetscOptionsEnd(); 65907a3e9cSStefano Zampini 66907a3e9cSStefano Zampini PetscCall(DMViewFromOptions(dm, NULL, "-dm_plex_dd_dm_view")); 67907a3e9cSStefano Zampini PetscCall(DMPlexDistributeOverlap_Internal(dm, ddovl + 1, PETSC_COMM_SELF, oname, &migrationSF, &odm)); 68907a3e9cSStefano Zampini if (!odm) PetscCall(DMClone(dm, &odm)); 69907a3e9cSStefano Zampini if (migrationSF) PetscCall(PetscSFViewFromOptions(migrationSF, (PetscObject)dm, "-dm_plex_dd_sf_view")); 70907a3e9cSStefano Zampini 71907a3e9cSStefano Zampini PetscCall(DMPlexGetMaxProjectionHeight(dm, &mh)); 72907a3e9cSStefano Zampini PetscCall(DMPlexSetMaxProjectionHeight(odm, mh)); 73907a3e9cSStefano Zampini 74907a3e9cSStefano Zampini /* share discretization */ 75907a3e9cSStefano Zampini /* TODO Labels for regions may need to updated, 76907a3e9cSStefano Zampini now it uses the original ones, not the ones for the odm. 77907a3e9cSStefano Zampini Not sure what to do here */ 78907a3e9cSStefano Zampini /* PetscCall(DMCopyDisc(dm, odm)); */ 79907a3e9cSStefano Zampini PetscCall(DMGetDS(odm, &ds)); 80907a3e9cSStefano Zampini if (!ds) { 81907a3e9cSStefano Zampini PetscCall(DMGetLocalSection(dm, &sec)); 82907a3e9cSStefano Zampini PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)odm), &tsec)); 83907a3e9cSStefano Zampini if (migrationSF) { 84907a3e9cSStefano Zampini PetscCall(PetscSFDistributeSection(migrationSF, sec, NULL, tsec)); 85907a3e9cSStefano Zampini } else { 86907a3e9cSStefano Zampini PetscCall(PetscSectionCopy(sec, tsec)); 87907a3e9cSStefano Zampini } 88907a3e9cSStefano Zampini PetscCall(DMSetLocalSection(dm, tsec)); 89907a3e9cSStefano Zampini PetscCall(PetscSectionDestroy(&tsec)); 90907a3e9cSStefano Zampini } 91e4d5475eSStefano Zampini /* TODO: what if these values changes? add to some DM hook? */ 92e4d5475eSStefano Zampini PetscCall(DMTransferMaterialParameters(dm, migrationSF, odm)); 93907a3e9cSStefano Zampini 94*0fc0a6c4SStefano Zampini PetscCall(DMViewFromOptions(odm, (PetscObject)dm, "-dm_plex_dd_overlap_dm_view")); 95907a3e9cSStefano Zampini #if 0 96907a3e9cSStefano Zampini { 97*0fc0a6c4SStefano Zampini DM seqdm; 98*0fc0a6c4SStefano Zampini Vec val; 99*0fc0a6c4SStefano Zampini IS is; 100*0fc0a6c4SStefano Zampini PetscInt vStart, vEnd; 101*0fc0a6c4SStefano Zampini const PetscInt *vnum; 102907a3e9cSStefano Zampini char name[256]; 103907a3e9cSStefano Zampini PetscViewer viewer; 104907a3e9cSStefano Zampini 105*0fc0a6c4SStefano Zampini PetscCall(DMPlexDistributeOverlap_Internal(dm, 0, PETSC_COMM_SELF, "local_mesh", NULL, &seqdm)); 106*0fc0a6c4SStefano Zampini PetscCall(PetscSNPrintf(name, sizeof(name), "local_mesh_%d.vtu", PetscGlobalRank)); 107*0fc0a6c4SStefano Zampini PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)seqdm), name, FILE_MODE_WRITE, &viewer)); 108*0fc0a6c4SStefano Zampini PetscCall(DMGetLabel(seqdm, "local_mesh", &label)); 109*0fc0a6c4SStefano Zampini PetscCall(DMPlexLabelComplete(seqdm, label)); 110*0fc0a6c4SStefano Zampini PetscCall(DMPlexCreateLabelField(seqdm, label, &val)); 111*0fc0a6c4SStefano Zampini PetscCall(VecView(val, viewer)); 112*0fc0a6c4SStefano Zampini PetscCall(VecDestroy(&val)); 113*0fc0a6c4SStefano Zampini PetscCall(PetscViewerDestroy(&viewer)); 114*0fc0a6c4SStefano Zampini 115*0fc0a6c4SStefano Zampini PetscCall(PetscSNPrintf(name, sizeof(name), "asm_vertices_%d.vtu", PetscGlobalRank)); 116*0fc0a6c4SStefano Zampini PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)seqdm), name, FILE_MODE_WRITE, &viewer)); 117*0fc0a6c4SStefano Zampini PetscCall(DMCreateLabel(seqdm, "asm_vertices")); 118*0fc0a6c4SStefano Zampini PetscCall(DMGetLabel(seqdm, "asm_vertices", &label)); 119*0fc0a6c4SStefano Zampini PetscCall(DMPlexGetVertexNumbering(dm, &is)); 120*0fc0a6c4SStefano Zampini PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 121*0fc0a6c4SStefano Zampini PetscCall(ISGetIndices(is, &vnum)); 122*0fc0a6c4SStefano Zampini for (PetscInt v = 0; v < vEnd - vStart; v++) { 123*0fc0a6c4SStefano Zampini if (vnum[v] < 0) continue; 124*0fc0a6c4SStefano Zampini PetscCall(DMLabelSetValue(label, v + vStart, 1)); 125*0fc0a6c4SStefano Zampini } 126*0fc0a6c4SStefano Zampini PetscCall(DMPlexCreateLabelField(seqdm, label, &val)); 127*0fc0a6c4SStefano Zampini PetscCall(VecView(val, viewer)); 128*0fc0a6c4SStefano Zampini PetscCall(VecDestroy(&val)); 129*0fc0a6c4SStefano Zampini PetscCall(ISRestoreIndices(is, &vnum)); 130*0fc0a6c4SStefano Zampini PetscCall(PetscViewerDestroy(&viewer)); 131*0fc0a6c4SStefano Zampini 132*0fc0a6c4SStefano Zampini PetscCall(DMDestroy(&seqdm)); 133907a3e9cSStefano Zampini PetscCall(PetscSNPrintf(name, sizeof(name), "ovl_mesh_%d.vtu", PetscGlobalRank)); 134907a3e9cSStefano Zampini PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)odm), name, FILE_MODE_WRITE, &viewer)); 135907a3e9cSStefano Zampini PetscCall(DMView(odm, viewer)); 136907a3e9cSStefano Zampini PetscCall(PetscViewerDestroy(&viewer)); 137907a3e9cSStefano Zampini } 138907a3e9cSStefano Zampini #endif 139907a3e9cSStefano Zampini 140907a3e9cSStefano Zampini /* propagate original global ordering to overlapping DM */ 141907a3e9cSStefano Zampini PetscCall(DMGetSectionSF(dm, §ionSF)); 142907a3e9cSStefano Zampini PetscCall(DMGetLocalSection(dm, &sec)); 143907a3e9cSStefano Zampini PetscCall(PetscSectionGetStorageSize(sec, &nl)); 144907a3e9cSStefano Zampini PetscCall(DMGetGlobalVector(dm, &gvec)); 145907a3e9cSStefano Zampini PetscCall(VecGetOwnershipRange(gvec, &rst, &ren)); 146907a3e9cSStefano Zampini PetscCall(DMRestoreGlobalVector(dm, &gvec)); 147907a3e9cSStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF, ren - rst, rst, 1, &gi_is)); /* non-overlapping dofs */ 148907a3e9cSStefano Zampini PetscCall(PetscMalloc1(nl, &lidxs)); 149907a3e9cSStefano Zampini for (PetscInt i = 0; i < nl; i++) lidxs[i] = -1; 150907a3e9cSStefano Zampini PetscCall(ISGetIndices(gi_is, (const PetscInt **)&gidxs)); 151907a3e9cSStefano Zampini PetscCall(PetscSFBcastBegin(sectionSF, MPIU_INT, gidxs, lidxs, MPI_REPLACE)); 152907a3e9cSStefano Zampini PetscCall(PetscSFBcastEnd(sectionSF, MPIU_INT, gidxs, lidxs, MPI_REPLACE)); 153907a3e9cSStefano Zampini PetscCall(ISRestoreIndices(gi_is, (const PetscInt **)&gidxs)); 154907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), nl, lidxs, PETSC_OWN_POINTER, &lis)); 155907a3e9cSStefano Zampini PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)odm), &tsec)); 156907a3e9cSStefano Zampini PetscCall(DMPlexDistributeFieldIS(dm, migrationSF, sec, lis, tsec, &gis)); 157907a3e9cSStefano Zampini PetscCall(PetscSectionDestroy(&tsec)); 158907a3e9cSStefano Zampini PetscCall(ISDestroy(&lis)); 159907a3e9cSStefano Zampini PetscCall(PetscSFDestroy(&migrationSF)); 160907a3e9cSStefano Zampini 161907a3e9cSStefano Zampini /* make dofs on the overlap boundary (not the global boundary) constrained */ 162907a3e9cSStefano Zampini PetscCall(DMGetLabel(odm, oname, &label)); 163907a3e9cSStefano Zampini PetscCall(DMPlexLabelComplete(odm, label)); 164907a3e9cSStefano Zampini PetscCall(DMGetLocalSection(odm, &tsec)); 165907a3e9cSStefano Zampini PetscCall(PetscSectionGetChart(tsec, &pStart, &pEnd)); 166907a3e9cSStefano Zampini PetscCall(DMLabelCreateIndex(label, pStart, pEnd)); 167907a3e9cSStefano Zampini PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)tsec), &sec)); 168907a3e9cSStefano Zampini PetscCall(PetscSectionCopy_Internal(tsec, sec, label->bt)); 169907a3e9cSStefano Zampini PetscCall(DMSetLocalSection(odm, sec)); 170907a3e9cSStefano Zampini PetscCall(PetscSectionDestroy(&sec)); 171907a3e9cSStefano Zampini PetscCall(DMRemoveLabel(odm, oname, NULL)); 172907a3e9cSStefano Zampini 173907a3e9cSStefano Zampini /* Create index sets for dofs in the overlap dm */ 174907a3e9cSStefano Zampini PetscCall(DMGetSectionSF(odm, §ionSF)); 175907a3e9cSStefano Zampini PetscCall(DMGetLocalSection(odm, &olsec)); 176907a3e9cSStefano Zampini PetscCall(DMGetGlobalSection(odm, &ogsec)); 177907a3e9cSStefano Zampini PetscCall(PetscSectionViewFromOptions(ogsec, (PetscObject)dm, "-dm_plex_dd_overlap_gsection_view")); 178907a3e9cSStefano Zampini PetscCall(PetscSectionViewFromOptions(olsec, (PetscObject)dm, "-dm_plex_dd_overlap_lsection_view")); 179907a3e9cSStefano Zampini ni = ren - rst; 180907a3e9cSStefano Zampini PetscCall(PetscSectionGetConstrainedStorageSize(ogsec, &no)); /* dofs in the overlap */ 181907a3e9cSStefano Zampini PetscCall(PetscSectionGetStorageSize(olsec, &nl)); /* local dofs in the overlap */ 182907a3e9cSStefano Zampini PetscCall(PetscMalloc1(no, &gidxs)); 183907a3e9cSStefano Zampini PetscCall(ISGetIndices(gis, (const PetscInt **)&lidxs)); 184907a3e9cSStefano Zampini PetscCall(PetscSFReduceBegin(sectionSF, MPIU_INT, lidxs, gidxs, MPI_REPLACE)); 185907a3e9cSStefano Zampini PetscCall(PetscSFReduceEnd(sectionSF, MPIU_INT, lidxs, gidxs, MPI_REPLACE)); 186907a3e9cSStefano Zampini PetscCall(ISRestoreIndices(gis, (const PetscInt **)&lidxs)); 187907a3e9cSStefano Zampini 188907a3e9cSStefano Zampini /* non-overlapping dofs */ 189907a3e9cSStefano Zampini PetscCall(PetscMalloc1(no, &lidxs)); 190907a3e9cSStefano Zampini c = 0; 191907a3e9cSStefano Zampini for (PetscInt i = 0; i < no; i++) { 192907a3e9cSStefano Zampini if (gidxs[i] >= rst && gidxs[i] < ren) lidxs[c++] = i; 193907a3e9cSStefano Zampini } 194907a3e9cSStefano Zampini PetscCheck(c == ni, PETSC_COMM_SELF, PETSC_ERR_PLIB, "%" PetscInt_FMT " != %" PetscInt_FMT "\n", c, ni); 195907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), ni, lidxs, PETSC_OWN_POINTER, &li_is)); 196907a3e9cSStefano Zampini 197907a3e9cSStefano Zampini /* global dofs in the overlap */ 198907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), no, gidxs, PETSC_OWN_POINTER, &go_is)); 199907a3e9cSStefano Zampini PetscCall(ISViewFromOptions(go_is, (PetscObject)dm, "-dm_plex_dd_overlap_gois_view")); 200907a3e9cSStefano Zampini /* PetscCall(ISCreateStride(PetscObjectComm((PetscObject)odm), no, 0, 1, &lo_is)); */ 201907a3e9cSStefano Zampini 202907a3e9cSStefano Zampini /* local dofs of the overlapping subdomain (we actually need only dofs on the boundary of the subdomain) */ 203907a3e9cSStefano Zampini PetscCall(PetscMalloc1(nl, &lidxs)); 204907a3e9cSStefano Zampini PetscCall(PetscMalloc1(nl, &gidxs)); 205907a3e9cSStefano Zampini PetscCall(ISGetIndices(gis, (const PetscInt **)&tidxs)); 206907a3e9cSStefano Zampini c = 0; 207907a3e9cSStefano Zampini for (PetscInt i = 0; i < nl; i++) { 208907a3e9cSStefano Zampini if (tidxs[i] < 0) continue; 209907a3e9cSStefano Zampini lidxs[c] = i; 210907a3e9cSStefano Zampini gidxs[c] = tidxs[i]; 211907a3e9cSStefano Zampini c++; 212907a3e9cSStefano Zampini } 213907a3e9cSStefano Zampini PetscCall(ISRestoreIndices(gis, (const PetscInt **)&tidxs)); 214907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), c, gidxs, PETSC_OWN_POINTER, &gl_is)); 215907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), c, lidxs, PETSC_OWN_POINTER, &ll_is)); 216907a3e9cSStefano Zampini PetscCall(ISViewFromOptions(gl_is, (PetscObject)dm, "-dm_plex_dd_overlap_glis_view")); 217907a3e9cSStefano Zampini 218907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_gi", (PetscObject)gi_is)); 219907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_li", (PetscObject)li_is)); 220907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_go", (PetscObject)go_is)); 221907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_gl", (PetscObject)gl_is)); 222907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_ll", (PetscObject)ll_is)); 223907a3e9cSStefano Zampini 224907a3e9cSStefano Zampini if (innerises) (*innerises)[0] = gi_is; 225907a3e9cSStefano Zampini else PetscCall(ISDestroy(&gi_is)); 226907a3e9cSStefano Zampini if (outerises) (*outerises)[0] = go_is; 227907a3e9cSStefano Zampini else PetscCall(ISDestroy(&go_is)); 228907a3e9cSStefano Zampini if (dms) (*dms)[0] = odm; 229907a3e9cSStefano Zampini else PetscCall(DMDestroy(&odm)); 230907a3e9cSStefano Zampini PetscCall(ISDestroy(&li_is)); 231907a3e9cSStefano Zampini PetscCall(ISDestroy(&gl_is)); 232907a3e9cSStefano Zampini PetscCall(ISDestroy(&ll_is)); 233907a3e9cSStefano Zampini PetscCall(ISDestroy(&gis)); 234907a3e9cSStefano Zampini 235907a3e9cSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 236907a3e9cSStefano Zampini } 237907a3e9cSStefano Zampini 238907a3e9cSStefano Zampini PetscErrorCode DMCreateDomainDecompositionScatters_Plex(DM dm, PetscInt n, DM *subdms, VecScatter **iscat, VecScatter **oscat, VecScatter **lscat) 239907a3e9cSStefano Zampini { 240907a3e9cSStefano Zampini Vec gvec, svec, lvec; 241e4d5475eSStefano Zampini IS gi_is, li_is, go_is, gl_is, ll_is; 242907a3e9cSStefano Zampini 243907a3e9cSStefano Zampini PetscFunctionBegin; 244907a3e9cSStefano Zampini if (iscat) PetscCall(PetscMalloc1(n, iscat)); 245907a3e9cSStefano Zampini if (oscat) PetscCall(PetscMalloc1(n, oscat)); 246907a3e9cSStefano Zampini if (lscat) PetscCall(PetscMalloc1(n, lscat)); 247907a3e9cSStefano Zampini 248907a3e9cSStefano Zampini PetscCall(DMGetGlobalVector(dm, &gvec)); 249907a3e9cSStefano Zampini for (PetscInt i = 0; i < n; i++) { 250907a3e9cSStefano Zampini PetscCall(DMGetGlobalVector(subdms[i], &svec)); 251907a3e9cSStefano Zampini PetscCall(DMGetLocalVector(subdms[i], &lvec)); 252907a3e9cSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_gi", (PetscObject *)&gi_is)); 253907a3e9cSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_li", (PetscObject *)&li_is)); 254907a3e9cSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_go", (PetscObject *)&go_is)); 255907a3e9cSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_gl", (PetscObject *)&gl_is)); 256907a3e9cSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_ll", (PetscObject *)&ll_is)); 257e4d5475eSStefano Zampini PetscCheck(gi_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition"); 258e4d5475eSStefano Zampini PetscCheck(li_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition"); 259e4d5475eSStefano Zampini PetscCheck(go_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition"); 260e4d5475eSStefano Zampini PetscCheck(gl_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition"); 261e4d5475eSStefano Zampini PetscCheck(ll_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition"); 262907a3e9cSStefano Zampini if (iscat) PetscCall(VecScatterCreate(gvec, gi_is, svec, li_is, &(*iscat)[i])); 263e4d5475eSStefano Zampini if (oscat) PetscCall(VecScatterCreate(gvec, go_is, svec, NULL, &(*oscat)[i])); 264907a3e9cSStefano Zampini if (lscat) PetscCall(VecScatterCreate(gvec, gl_is, lvec, ll_is, &(*lscat)[i])); 265907a3e9cSStefano Zampini PetscCall(DMRestoreGlobalVector(subdms[i], &svec)); 266907a3e9cSStefano Zampini PetscCall(DMRestoreLocalVector(subdms[i], &lvec)); 267907a3e9cSStefano Zampini } 268907a3e9cSStefano Zampini PetscCall(DMRestoreGlobalVector(dm, &gvec)); 269907a3e9cSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 270907a3e9cSStefano Zampini } 271907a3e9cSStefano Zampini 272907a3e9cSStefano Zampini /* 273907a3e9cSStefano Zampini DMCreateNeumannOverlap_Plex - Generates an IS, an unassembled (Neumann) Mat, a setup function, and the corresponding context to be used by PCHPDDM. 274907a3e9cSStefano Zampini 275907a3e9cSStefano Zampini Input Parameter: 276907a3e9cSStefano Zampini . dm - preconditioner context 277907a3e9cSStefano Zampini 278907a3e9cSStefano Zampini Output Parameters: 279907a3e9cSStefano Zampini + ovl - index set of overlapping subdomains 280907a3e9cSStefano Zampini . J - unassembled (Neumann) local matrix 281907a3e9cSStefano Zampini . setup - function for generating the matrix 282907a3e9cSStefano Zampini - setup_ctx - context for setup 283907a3e9cSStefano Zampini 284907a3e9cSStefano Zampini Options Database Keys: 285907a3e9cSStefano Zampini + -dm_plex_view_neumann_original - view the DM without overlap 286907a3e9cSStefano Zampini - -dm_plex_view_neumann_overlap - view the DM with overlap as needed by PCHPDDM 287907a3e9cSStefano Zampini 288907a3e9cSStefano Zampini Level: advanced 289907a3e9cSStefano Zampini 290907a3e9cSStefano Zampini .seealso: `DMCreate()`, `DM`, `MATIS`, `PCHPDDM`, `PCHPDDMSetAuxiliaryMat()` 291907a3e9cSStefano Zampini */ 292907a3e9cSStefano Zampini PetscErrorCode DMCreateNeumannOverlap_Plex(DM dm, IS *ovl, Mat *J, PetscErrorCode (**setup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void **setup_ctx) 293907a3e9cSStefano Zampini { 294907a3e9cSStefano Zampini DM odm; 295907a3e9cSStefano Zampini Mat pJ; 296907a3e9cSStefano Zampini PetscSF sf = NULL; 297907a3e9cSStefano Zampini PetscSection sec, osec; 298907a3e9cSStefano Zampini ISLocalToGlobalMapping l2g; 299907a3e9cSStefano Zampini const PetscInt *idxs; 300907a3e9cSStefano Zampini PetscInt n, mh; 301907a3e9cSStefano Zampini 302907a3e9cSStefano Zampini PetscFunctionBegin; 303907a3e9cSStefano Zampini *setup = NULL; 304907a3e9cSStefano Zampini *setup_ctx = NULL; 305907a3e9cSStefano Zampini *ovl = NULL; 306907a3e9cSStefano Zampini *J = NULL; 307907a3e9cSStefano Zampini 308907a3e9cSStefano Zampini /* Overlapped mesh 309907a3e9cSStefano Zampini overlap is a little more generous, since it is not computed starting from the owned (Dirichlet) points, but from the locally owned cells */ 310907a3e9cSStefano Zampini PetscCall(DMPlexDistributeOverlap(dm, 1, &sf, &odm)); 311907a3e9cSStefano Zampini if (!odm) { 312907a3e9cSStefano Zampini PetscCall(PetscSFDestroy(&sf)); 313907a3e9cSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 314907a3e9cSStefano Zampini } 315907a3e9cSStefano Zampini 316907a3e9cSStefano Zampini /* share discretization */ 317907a3e9cSStefano Zampini PetscCall(DMGetLocalSection(dm, &sec)); 318907a3e9cSStefano Zampini PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)sec), &osec)); 319907a3e9cSStefano Zampini PetscCall(PetscSFDistributeSection(sf, sec, NULL, osec)); 320907a3e9cSStefano Zampini /* what to do here? using both is fine? */ 321907a3e9cSStefano Zampini PetscCall(DMSetLocalSection(odm, osec)); 322907a3e9cSStefano Zampini PetscCall(DMCopyDisc(dm, odm)); 323907a3e9cSStefano Zampini PetscCall(DMPlexGetMaxProjectionHeight(dm, &mh)); 324907a3e9cSStefano Zampini PetscCall(DMPlexSetMaxProjectionHeight(odm, mh)); 325907a3e9cSStefano Zampini PetscCall(PetscSectionDestroy(&osec)); 326907a3e9cSStefano Zampini 327907a3e9cSStefano Zampini /* material parameters */ 328e4d5475eSStefano Zampini /* TODO: what if these values changes? add to some DM hook? */ 329e4d5475eSStefano Zampini PetscCall(DMTransferMaterialParameters(dm, sf, odm)); 330907a3e9cSStefano Zampini PetscCall(PetscSFDestroy(&sf)); 331907a3e9cSStefano Zampini 332907a3e9cSStefano Zampini PetscCall(DMViewFromOptions(dm, NULL, "-dm_plex_view_neumann_original")); 333907a3e9cSStefano Zampini PetscCall(PetscObjectSetName((PetscObject)odm, "OVL")); 334907a3e9cSStefano Zampini PetscCall(DMViewFromOptions(odm, NULL, "-dm_plex_view_neumann_overlap")); 335907a3e9cSStefano Zampini 336907a3e9cSStefano Zampini /* MATIS for the overlap region 337907a3e9cSStefano Zampini the HPDDM interface wants local matrices, 338907a3e9cSStefano Zampini we attach the global MATIS to the overlap IS, 339907a3e9cSStefano Zampini since we need it to do assembly */ 340907a3e9cSStefano Zampini PetscCall(DMSetMatType(odm, MATIS)); 341907a3e9cSStefano Zampini PetscCall(DMCreateMatrix(odm, &pJ)); 342907a3e9cSStefano Zampini PetscCall(MatISGetLocalMat(pJ, J)); 343907a3e9cSStefano Zampini PetscCall(PetscObjectReference((PetscObject)*J)); 344907a3e9cSStefano Zampini 345907a3e9cSStefano Zampini /* overlap IS */ 346907a3e9cSStefano Zampini PetscCall(MatISGetLocalToGlobalMapping(pJ, &l2g, NULL)); 347907a3e9cSStefano Zampini PetscCall(ISLocalToGlobalMappingGetSize(l2g, &n)); 348907a3e9cSStefano Zampini PetscCall(ISLocalToGlobalMappingGetIndices(l2g, &idxs)); 349907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), n, idxs, PETSC_COPY_VALUES, ovl)); 350907a3e9cSStefano Zampini PetscCall(ISLocalToGlobalMappingRestoreIndices(l2g, &idxs)); 351907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)*ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject)pJ)); 352907a3e9cSStefano Zampini PetscCall(DMDestroy(&odm)); 353907a3e9cSStefano Zampini PetscCall(MatDestroy(&pJ)); 354907a3e9cSStefano Zampini 355907a3e9cSStefano Zampini /* special purpose setup function (composed in DMPlexSetSNESLocalFEM) */ 356907a3e9cSStefano Zampini PetscCall(PetscObjectQueryFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", setup)); 357907a3e9cSStefano Zampini if (*setup) PetscCall(PetscObjectCompose((PetscObject)*ovl, "_DM_Original_HPDDM", (PetscObject)dm)); 358907a3e9cSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 359907a3e9cSStefano Zampini } 360