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
DMTransferMaterialParameters(DM dm,PetscSF sf,DM odm)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));
25*a0d6ca37SStefano Zampini if (sf) {
26e4d5475eSStefano Zampini PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)sec), &osec));
27e4d5475eSStefano Zampini PetscCall(VecCreate(PetscObjectComm((PetscObject)A), &oAt));
28e4d5475eSStefano Zampini PetscCall(DMPlexDistributeField(dmAux, sf, sec, A, osec, oAt));
29*a0d6ca37SStefano Zampini } else {
30*a0d6ca37SStefano Zampini PetscCall(PetscObjectReference((PetscObject)sec));
31*a0d6ca37SStefano Zampini osec = sec;
32*a0d6ca37SStefano Zampini PetscCall(PetscObjectReference((PetscObject)A));
33*a0d6ca37SStefano Zampini oAt = A;
34*a0d6ca37SStefano Zampini }
35e4d5475eSStefano Zampini PetscCall(DMSetLocalSection(odmAux, osec));
36e4d5475eSStefano Zampini PetscCall(PetscSectionDestroy(&osec));
37e4d5475eSStefano Zampini PetscCall(DMCreateLocalVector(odmAux, &oA));
38e4d5475eSStefano Zampini PetscCall(DMDestroy(&odmAux));
39e4d5475eSStefano Zampini PetscCall(VecCopy(oAt, oA));
40e4d5475eSStefano Zampini PetscCall(VecDestroy(&oAt));
41e4d5475eSStefano Zampini PetscCall(DMSetAuxiliaryVec(odm, NULL, 0, 0, oA));
42e4d5475eSStefano Zampini PetscCall(VecDestroy(&oA));
43e4d5475eSStefano Zampini }
44e4d5475eSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
45e4d5475eSStefano Zampini }
46e4d5475eSStefano Zampini
DMCreateDomainDecomposition_Plex(DM dm,PetscInt * nsub,char *** names,IS ** innerises,IS ** outerises,DM ** dms)47907a3e9cSStefano Zampini PetscErrorCode DMCreateDomainDecomposition_Plex(DM dm, PetscInt *nsub, char ***names, IS **innerises, IS **outerises, DM **dms)
48907a3e9cSStefano Zampini {
49907a3e9cSStefano Zampini DM odm;
50907a3e9cSStefano Zampini PetscSF migrationSF = NULL, sectionSF;
51907a3e9cSStefano Zampini PetscSection sec, tsec, ogsec, olsec;
52907a3e9cSStefano Zampini PetscInt n, mh, ddovl = 0, pStart, pEnd, ni, no, nl;
53907a3e9cSStefano Zampini PetscDS ds;
54907a3e9cSStefano Zampini DMLabel label;
55907a3e9cSStefano Zampini const char *oname = "__internal_plex_dd_ovl_";
56907a3e9cSStefano Zampini IS gi_is, li_is, go_is, gl_is, ll_is;
57907a3e9cSStefano Zampini IS gis, lis;
58907a3e9cSStefano Zampini PetscInt rst, ren, c, *gidxs, *lidxs, *tidxs;
59907a3e9cSStefano Zampini Vec gvec;
60907a3e9cSStefano Zampini
61907a3e9cSStefano Zampini PetscFunctionBegin;
62907a3e9cSStefano Zampini n = 1;
63907a3e9cSStefano Zampini if (nsub) *nsub = n;
64907a3e9cSStefano Zampini if (names) PetscCall(PetscCalloc1(n, names));
65907a3e9cSStefano Zampini if (innerises) PetscCall(PetscCalloc1(n, innerises));
66907a3e9cSStefano Zampini if (outerises) PetscCall(PetscCalloc1(n, outerises));
67907a3e9cSStefano Zampini if (dms) PetscCall(PetscCalloc1(n, dms));
68907a3e9cSStefano Zampini
69907a3e9cSStefano Zampini PetscObjectOptionsBegin((PetscObject)dm);
70907a3e9cSStefano Zampini PetscCall(PetscOptionsBoundedInt("-dm_plex_dd_overlap", "The size of the overlap halo for the subdomains", "DMCreateDomainDecomposition", ddovl, &ddovl, NULL, 0));
71907a3e9cSStefano Zampini PetscOptionsEnd();
72907a3e9cSStefano Zampini
73907a3e9cSStefano Zampini PetscCall(DMViewFromOptions(dm, NULL, "-dm_plex_dd_dm_view"));
74907a3e9cSStefano Zampini PetscCall(DMPlexDistributeOverlap_Internal(dm, ddovl + 1, PETSC_COMM_SELF, oname, &migrationSF, &odm));
75907a3e9cSStefano Zampini if (!odm) PetscCall(DMClone(dm, &odm));
76907a3e9cSStefano Zampini if (migrationSF) PetscCall(PetscSFViewFromOptions(migrationSF, (PetscObject)dm, "-dm_plex_dd_sf_view"));
77907a3e9cSStefano Zampini
78907a3e9cSStefano Zampini PetscCall(DMPlexGetMaxProjectionHeight(dm, &mh));
79907a3e9cSStefano Zampini PetscCall(DMPlexSetMaxProjectionHeight(odm, mh));
80907a3e9cSStefano Zampini
81907a3e9cSStefano Zampini /* share discretization */
82907a3e9cSStefano Zampini /* TODO Labels for regions may need to updated,
83907a3e9cSStefano Zampini now it uses the original ones, not the ones for the odm.
84907a3e9cSStefano Zampini Not sure what to do here */
85907a3e9cSStefano Zampini /* PetscCall(DMCopyDisc(dm, odm)); */
86907a3e9cSStefano Zampini PetscCall(DMGetDS(odm, &ds));
87907a3e9cSStefano Zampini if (!ds) {
88907a3e9cSStefano Zampini PetscCall(DMGetLocalSection(dm, &sec));
89907a3e9cSStefano Zampini PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)odm), &tsec));
90907a3e9cSStefano Zampini if (migrationSF) {
91907a3e9cSStefano Zampini PetscCall(PetscSFDistributeSection(migrationSF, sec, NULL, tsec));
92907a3e9cSStefano Zampini } else {
93907a3e9cSStefano Zampini PetscCall(PetscSectionCopy(sec, tsec));
94907a3e9cSStefano Zampini }
95907a3e9cSStefano Zampini PetscCall(DMSetLocalSection(dm, tsec));
96907a3e9cSStefano Zampini PetscCall(PetscSectionDestroy(&tsec));
97907a3e9cSStefano Zampini }
98e4d5475eSStefano Zampini /* TODO: what if these values changes? add to some DM hook? */
99e4d5475eSStefano Zampini PetscCall(DMTransferMaterialParameters(dm, migrationSF, odm));
100907a3e9cSStefano Zampini
1010fc0a6c4SStefano Zampini PetscCall(DMViewFromOptions(odm, (PetscObject)dm, "-dm_plex_dd_overlap_dm_view"));
102907a3e9cSStefano Zampini #if 0
103907a3e9cSStefano Zampini {
1040fc0a6c4SStefano Zampini DM seqdm;
1050fc0a6c4SStefano Zampini Vec val;
1060fc0a6c4SStefano Zampini IS is;
1070fc0a6c4SStefano Zampini PetscInt vStart, vEnd;
1080fc0a6c4SStefano Zampini const PetscInt *vnum;
109907a3e9cSStefano Zampini char name[256];
110907a3e9cSStefano Zampini PetscViewer viewer;
111907a3e9cSStefano Zampini
1120fc0a6c4SStefano Zampini PetscCall(DMPlexDistributeOverlap_Internal(dm, 0, PETSC_COMM_SELF, "local_mesh", NULL, &seqdm));
1130fc0a6c4SStefano Zampini PetscCall(PetscSNPrintf(name, sizeof(name), "local_mesh_%d.vtu", PetscGlobalRank));
1140fc0a6c4SStefano Zampini PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)seqdm), name, FILE_MODE_WRITE, &viewer));
1150fc0a6c4SStefano Zampini PetscCall(DMGetLabel(seqdm, "local_mesh", &label));
1160fc0a6c4SStefano Zampini PetscCall(DMPlexLabelComplete(seqdm, label));
1170fc0a6c4SStefano Zampini PetscCall(DMPlexCreateLabelField(seqdm, label, &val));
1180fc0a6c4SStefano Zampini PetscCall(VecView(val, viewer));
1190fc0a6c4SStefano Zampini PetscCall(VecDestroy(&val));
1200fc0a6c4SStefano Zampini PetscCall(PetscViewerDestroy(&viewer));
1210fc0a6c4SStefano Zampini
1220fc0a6c4SStefano Zampini PetscCall(PetscSNPrintf(name, sizeof(name), "asm_vertices_%d.vtu", PetscGlobalRank));
1230fc0a6c4SStefano Zampini PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)seqdm), name, FILE_MODE_WRITE, &viewer));
1240fc0a6c4SStefano Zampini PetscCall(DMCreateLabel(seqdm, "asm_vertices"));
1250fc0a6c4SStefano Zampini PetscCall(DMGetLabel(seqdm, "asm_vertices", &label));
1260fc0a6c4SStefano Zampini PetscCall(DMPlexGetVertexNumbering(dm, &is));
1270fc0a6c4SStefano Zampini PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1280fc0a6c4SStefano Zampini PetscCall(ISGetIndices(is, &vnum));
1290fc0a6c4SStefano Zampini for (PetscInt v = 0; v < vEnd - vStart; v++) {
1300fc0a6c4SStefano Zampini if (vnum[v] < 0) continue;
1310fc0a6c4SStefano Zampini PetscCall(DMLabelSetValue(label, v + vStart, 1));
1320fc0a6c4SStefano Zampini }
1330fc0a6c4SStefano Zampini PetscCall(DMPlexCreateLabelField(seqdm, label, &val));
1340fc0a6c4SStefano Zampini PetscCall(VecView(val, viewer));
1350fc0a6c4SStefano Zampini PetscCall(VecDestroy(&val));
1360fc0a6c4SStefano Zampini PetscCall(ISRestoreIndices(is, &vnum));
1370fc0a6c4SStefano Zampini PetscCall(PetscViewerDestroy(&viewer));
1380fc0a6c4SStefano Zampini
1390fc0a6c4SStefano Zampini PetscCall(DMDestroy(&seqdm));
140907a3e9cSStefano Zampini PetscCall(PetscSNPrintf(name, sizeof(name), "ovl_mesh_%d.vtu", PetscGlobalRank));
141907a3e9cSStefano Zampini PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)odm), name, FILE_MODE_WRITE, &viewer));
142907a3e9cSStefano Zampini PetscCall(DMView(odm, viewer));
143907a3e9cSStefano Zampini PetscCall(PetscViewerDestroy(&viewer));
144907a3e9cSStefano Zampini }
145907a3e9cSStefano Zampini #endif
146907a3e9cSStefano Zampini
147907a3e9cSStefano Zampini /* propagate original global ordering to overlapping DM */
148907a3e9cSStefano Zampini PetscCall(DMGetSectionSF(dm, §ionSF));
149907a3e9cSStefano Zampini PetscCall(DMGetLocalSection(dm, &sec));
150907a3e9cSStefano Zampini PetscCall(PetscSectionGetStorageSize(sec, &nl));
151907a3e9cSStefano Zampini PetscCall(DMGetGlobalVector(dm, &gvec));
152907a3e9cSStefano Zampini PetscCall(VecGetOwnershipRange(gvec, &rst, &ren));
153907a3e9cSStefano Zampini PetscCall(DMRestoreGlobalVector(dm, &gvec));
154907a3e9cSStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF, ren - rst, rst, 1, &gi_is)); /* non-overlapping dofs */
155907a3e9cSStefano Zampini PetscCall(PetscMalloc1(nl, &lidxs));
156907a3e9cSStefano Zampini for (PetscInt i = 0; i < nl; i++) lidxs[i] = -1;
157907a3e9cSStefano Zampini PetscCall(ISGetIndices(gi_is, (const PetscInt **)&gidxs));
158907a3e9cSStefano Zampini PetscCall(PetscSFBcastBegin(sectionSF, MPIU_INT, gidxs, lidxs, MPI_REPLACE));
159907a3e9cSStefano Zampini PetscCall(PetscSFBcastEnd(sectionSF, MPIU_INT, gidxs, lidxs, MPI_REPLACE));
160907a3e9cSStefano Zampini PetscCall(ISRestoreIndices(gi_is, (const PetscInt **)&gidxs));
161907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), nl, lidxs, PETSC_OWN_POINTER, &lis));
162*a0d6ca37SStefano Zampini if (migrationSF) {
163907a3e9cSStefano Zampini PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)odm), &tsec));
164907a3e9cSStefano Zampini PetscCall(DMPlexDistributeFieldIS(dm, migrationSF, sec, lis, tsec, &gis));
165*a0d6ca37SStefano Zampini } else {
166*a0d6ca37SStefano Zampini PetscCall(PetscObjectReference((PetscObject)lis));
167*a0d6ca37SStefano Zampini gis = lis;
168*a0d6ca37SStefano Zampini tsec = NULL;
169*a0d6ca37SStefano Zampini }
170907a3e9cSStefano Zampini PetscCall(PetscSectionDestroy(&tsec));
171907a3e9cSStefano Zampini PetscCall(ISDestroy(&lis));
172907a3e9cSStefano Zampini PetscCall(PetscSFDestroy(&migrationSF));
173907a3e9cSStefano Zampini
174907a3e9cSStefano Zampini /* make dofs on the overlap boundary (not the global boundary) constrained */
175907a3e9cSStefano Zampini PetscCall(DMGetLabel(odm, oname, &label));
176*a0d6ca37SStefano Zampini if (label) {
177907a3e9cSStefano Zampini PetscCall(DMPlexLabelComplete(odm, label));
178907a3e9cSStefano Zampini PetscCall(DMGetLocalSection(odm, &tsec));
179907a3e9cSStefano Zampini PetscCall(PetscSectionGetChart(tsec, &pStart, &pEnd));
180907a3e9cSStefano Zampini PetscCall(DMLabelCreateIndex(label, pStart, pEnd));
181907a3e9cSStefano Zampini PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)tsec), &sec));
182907a3e9cSStefano Zampini PetscCall(PetscSectionCopy_Internal(tsec, sec, label->bt));
183907a3e9cSStefano Zampini PetscCall(DMSetLocalSection(odm, sec));
184907a3e9cSStefano Zampini PetscCall(PetscSectionDestroy(&sec));
185907a3e9cSStefano Zampini PetscCall(DMRemoveLabel(odm, oname, NULL));
186*a0d6ca37SStefano Zampini } else { /* sequential case */
187*a0d6ca37SStefano Zampini PetscCall(DMGetLocalSection(dm, &sec));
188*a0d6ca37SStefano Zampini PetscCall(DMSetLocalSection(odm, sec));
189*a0d6ca37SStefano Zampini }
190907a3e9cSStefano Zampini
191907a3e9cSStefano Zampini /* Create index sets for dofs in the overlap dm */
192907a3e9cSStefano Zampini PetscCall(DMGetSectionSF(odm, §ionSF));
193907a3e9cSStefano Zampini PetscCall(DMGetLocalSection(odm, &olsec));
194907a3e9cSStefano Zampini PetscCall(DMGetGlobalSection(odm, &ogsec));
195907a3e9cSStefano Zampini PetscCall(PetscSectionViewFromOptions(ogsec, (PetscObject)dm, "-dm_plex_dd_overlap_gsection_view"));
196907a3e9cSStefano Zampini PetscCall(PetscSectionViewFromOptions(olsec, (PetscObject)dm, "-dm_plex_dd_overlap_lsection_view"));
197907a3e9cSStefano Zampini ni = ren - rst;
198907a3e9cSStefano Zampini PetscCall(PetscSectionGetConstrainedStorageSize(ogsec, &no)); /* dofs in the overlap */
199907a3e9cSStefano Zampini PetscCall(PetscSectionGetStorageSize(olsec, &nl)); /* local dofs in the overlap */
200907a3e9cSStefano Zampini PetscCall(PetscMalloc1(no, &gidxs));
201907a3e9cSStefano Zampini PetscCall(ISGetIndices(gis, (const PetscInt **)&lidxs));
202907a3e9cSStefano Zampini PetscCall(PetscSFReduceBegin(sectionSF, MPIU_INT, lidxs, gidxs, MPI_REPLACE));
203907a3e9cSStefano Zampini PetscCall(PetscSFReduceEnd(sectionSF, MPIU_INT, lidxs, gidxs, MPI_REPLACE));
204907a3e9cSStefano Zampini PetscCall(ISRestoreIndices(gis, (const PetscInt **)&lidxs));
205907a3e9cSStefano Zampini
206907a3e9cSStefano Zampini /* non-overlapping dofs */
207907a3e9cSStefano Zampini PetscCall(PetscMalloc1(no, &lidxs));
208907a3e9cSStefano Zampini c = 0;
209907a3e9cSStefano Zampini for (PetscInt i = 0; i < no; i++) {
210907a3e9cSStefano Zampini if (gidxs[i] >= rst && gidxs[i] < ren) lidxs[c++] = i;
211907a3e9cSStefano Zampini }
212*a0d6ca37SStefano Zampini PetscCheck(c == ni, PETSC_COMM_SELF, PETSC_ERR_PLIB, "%" PetscInt_FMT " != %" PetscInt_FMT, c, ni);
213907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), ni, lidxs, PETSC_OWN_POINTER, &li_is));
214907a3e9cSStefano Zampini
215907a3e9cSStefano Zampini /* global dofs in the overlap */
216907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), no, gidxs, PETSC_OWN_POINTER, &go_is));
217907a3e9cSStefano Zampini PetscCall(ISViewFromOptions(go_is, (PetscObject)dm, "-dm_plex_dd_overlap_gois_view"));
218907a3e9cSStefano Zampini /* PetscCall(ISCreateStride(PetscObjectComm((PetscObject)odm), no, 0, 1, &lo_is)); */
219907a3e9cSStefano Zampini
220907a3e9cSStefano Zampini /* local dofs of the overlapping subdomain (we actually need only dofs on the boundary of the subdomain) */
221907a3e9cSStefano Zampini PetscCall(PetscMalloc1(nl, &lidxs));
222907a3e9cSStefano Zampini PetscCall(PetscMalloc1(nl, &gidxs));
223907a3e9cSStefano Zampini PetscCall(ISGetIndices(gis, (const PetscInt **)&tidxs));
224907a3e9cSStefano Zampini c = 0;
225907a3e9cSStefano Zampini for (PetscInt i = 0; i < nl; i++) {
226907a3e9cSStefano Zampini if (tidxs[i] < 0) continue;
227907a3e9cSStefano Zampini lidxs[c] = i;
228907a3e9cSStefano Zampini gidxs[c] = tidxs[i];
229907a3e9cSStefano Zampini c++;
230907a3e9cSStefano Zampini }
231907a3e9cSStefano Zampini PetscCall(ISRestoreIndices(gis, (const PetscInt **)&tidxs));
232907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), c, gidxs, PETSC_OWN_POINTER, &gl_is));
233907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), c, lidxs, PETSC_OWN_POINTER, &ll_is));
234907a3e9cSStefano Zampini PetscCall(ISViewFromOptions(gl_is, (PetscObject)dm, "-dm_plex_dd_overlap_glis_view"));
235907a3e9cSStefano Zampini
236907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_gi", (PetscObject)gi_is));
237907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_li", (PetscObject)li_is));
238907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_go", (PetscObject)go_is));
239907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_gl", (PetscObject)gl_is));
240907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_ll", (PetscObject)ll_is));
241907a3e9cSStefano Zampini
242907a3e9cSStefano Zampini if (innerises) (*innerises)[0] = gi_is;
243907a3e9cSStefano Zampini else PetscCall(ISDestroy(&gi_is));
244907a3e9cSStefano Zampini if (outerises) (*outerises)[0] = go_is;
245907a3e9cSStefano Zampini else PetscCall(ISDestroy(&go_is));
246907a3e9cSStefano Zampini if (dms) (*dms)[0] = odm;
247907a3e9cSStefano Zampini else PetscCall(DMDestroy(&odm));
248907a3e9cSStefano Zampini PetscCall(ISDestroy(&li_is));
249907a3e9cSStefano Zampini PetscCall(ISDestroy(&gl_is));
250907a3e9cSStefano Zampini PetscCall(ISDestroy(&ll_is));
251907a3e9cSStefano Zampini PetscCall(ISDestroy(&gis));
252907a3e9cSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
253907a3e9cSStefano Zampini }
254907a3e9cSStefano Zampini
DMCreateDomainDecompositionScatters_Plex(DM dm,PetscInt n,DM * subdms,VecScatter ** iscat,VecScatter ** oscat,VecScatter ** lscat)255907a3e9cSStefano Zampini PetscErrorCode DMCreateDomainDecompositionScatters_Plex(DM dm, PetscInt n, DM *subdms, VecScatter **iscat, VecScatter **oscat, VecScatter **lscat)
256907a3e9cSStefano Zampini {
257907a3e9cSStefano Zampini Vec gvec, svec, lvec;
258e4d5475eSStefano Zampini IS gi_is, li_is, go_is, gl_is, ll_is;
259907a3e9cSStefano Zampini
260907a3e9cSStefano Zampini PetscFunctionBegin;
261907a3e9cSStefano Zampini if (iscat) PetscCall(PetscMalloc1(n, iscat));
262907a3e9cSStefano Zampini if (oscat) PetscCall(PetscMalloc1(n, oscat));
263907a3e9cSStefano Zampini if (lscat) PetscCall(PetscMalloc1(n, lscat));
264907a3e9cSStefano Zampini
265907a3e9cSStefano Zampini PetscCall(DMGetGlobalVector(dm, &gvec));
266907a3e9cSStefano Zampini for (PetscInt i = 0; i < n; i++) {
267907a3e9cSStefano Zampini PetscCall(DMGetGlobalVector(subdms[i], &svec));
268907a3e9cSStefano Zampini PetscCall(DMGetLocalVector(subdms[i], &lvec));
269907a3e9cSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_gi", (PetscObject *)&gi_is));
270907a3e9cSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_li", (PetscObject *)&li_is));
271907a3e9cSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_go", (PetscObject *)&go_is));
272907a3e9cSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_gl", (PetscObject *)&gl_is));
273907a3e9cSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_ll", (PetscObject *)&ll_is));
274e4d5475eSStefano Zampini PetscCheck(gi_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition");
275e4d5475eSStefano Zampini PetscCheck(li_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition");
276e4d5475eSStefano Zampini PetscCheck(go_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition");
277e4d5475eSStefano Zampini PetscCheck(gl_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition");
278e4d5475eSStefano Zampini PetscCheck(ll_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition");
279907a3e9cSStefano Zampini if (iscat) PetscCall(VecScatterCreate(gvec, gi_is, svec, li_is, &(*iscat)[i]));
280e4d5475eSStefano Zampini if (oscat) PetscCall(VecScatterCreate(gvec, go_is, svec, NULL, &(*oscat)[i]));
281907a3e9cSStefano Zampini if (lscat) PetscCall(VecScatterCreate(gvec, gl_is, lvec, ll_is, &(*lscat)[i]));
282907a3e9cSStefano Zampini PetscCall(DMRestoreGlobalVector(subdms[i], &svec));
283907a3e9cSStefano Zampini PetscCall(DMRestoreLocalVector(subdms[i], &lvec));
284907a3e9cSStefano Zampini }
285907a3e9cSStefano Zampini PetscCall(DMRestoreGlobalVector(dm, &gvec));
286907a3e9cSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
287907a3e9cSStefano Zampini }
288907a3e9cSStefano Zampini
289907a3e9cSStefano Zampini /*
290907a3e9cSStefano Zampini DMCreateNeumannOverlap_Plex - Generates an IS, an unassembled (Neumann) Mat, a setup function, and the corresponding context to be used by PCHPDDM.
291907a3e9cSStefano Zampini
292907a3e9cSStefano Zampini Input Parameter:
293907a3e9cSStefano Zampini . dm - preconditioner context
294907a3e9cSStefano Zampini
295907a3e9cSStefano Zampini Output Parameters:
296907a3e9cSStefano Zampini + ovl - index set of overlapping subdomains
297907a3e9cSStefano Zampini . J - unassembled (Neumann) local matrix
298907a3e9cSStefano Zampini . setup - function for generating the matrix
299907a3e9cSStefano Zampini - setup_ctx - context for setup
300907a3e9cSStefano Zampini
301907a3e9cSStefano Zampini Options Database Keys:
302907a3e9cSStefano Zampini + -dm_plex_view_neumann_original - view the DM without overlap
303907a3e9cSStefano Zampini - -dm_plex_view_neumann_overlap - view the DM with overlap as needed by PCHPDDM
304907a3e9cSStefano Zampini
305907a3e9cSStefano Zampini Level: advanced
306907a3e9cSStefano Zampini
307907a3e9cSStefano Zampini .seealso: `DMCreate()`, `DM`, `MATIS`, `PCHPDDM`, `PCHPDDMSetAuxiliaryMat()`
308907a3e9cSStefano Zampini */
DMCreateNeumannOverlap_Plex(DM dm,IS * ovl,Mat * J,PetscErrorCode (** setup)(Mat,PetscReal,Vec,Vec,PetscReal,IS,void *),void ** setup_ctx)309907a3e9cSStefano Zampini PetscErrorCode DMCreateNeumannOverlap_Plex(DM dm, IS *ovl, Mat *J, PetscErrorCode (**setup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void **setup_ctx)
310907a3e9cSStefano Zampini {
311907a3e9cSStefano Zampini DM odm;
312907a3e9cSStefano Zampini Mat pJ;
313907a3e9cSStefano Zampini PetscSF sf = NULL;
314907a3e9cSStefano Zampini PetscSection sec, osec;
315907a3e9cSStefano Zampini ISLocalToGlobalMapping l2g;
316907a3e9cSStefano Zampini const PetscInt *idxs;
317907a3e9cSStefano Zampini PetscInt n, mh;
318907a3e9cSStefano Zampini
319907a3e9cSStefano Zampini PetscFunctionBegin;
320907a3e9cSStefano Zampini *setup = NULL;
321907a3e9cSStefano Zampini *setup_ctx = NULL;
322907a3e9cSStefano Zampini *ovl = NULL;
323907a3e9cSStefano Zampini *J = NULL;
324907a3e9cSStefano Zampini
325907a3e9cSStefano Zampini /* Overlapped mesh
326907a3e9cSStefano Zampini overlap is a little more generous, since it is not computed starting from the owned (Dirichlet) points, but from the locally owned cells */
327907a3e9cSStefano Zampini PetscCall(DMPlexDistributeOverlap(dm, 1, &sf, &odm));
328907a3e9cSStefano Zampini if (!odm) {
329907a3e9cSStefano Zampini PetscCall(PetscSFDestroy(&sf));
330907a3e9cSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
331907a3e9cSStefano Zampini }
332907a3e9cSStefano Zampini
333907a3e9cSStefano Zampini /* share discretization */
334907a3e9cSStefano Zampini PetscCall(DMGetLocalSection(dm, &sec));
335907a3e9cSStefano Zampini PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)sec), &osec));
336907a3e9cSStefano Zampini PetscCall(PetscSFDistributeSection(sf, sec, NULL, osec));
337907a3e9cSStefano Zampini /* what to do here? using both is fine? */
338907a3e9cSStefano Zampini PetscCall(DMSetLocalSection(odm, osec));
339907a3e9cSStefano Zampini PetscCall(DMCopyDisc(dm, odm));
340907a3e9cSStefano Zampini PetscCall(DMPlexGetMaxProjectionHeight(dm, &mh));
341907a3e9cSStefano Zampini PetscCall(DMPlexSetMaxProjectionHeight(odm, mh));
342907a3e9cSStefano Zampini PetscCall(PetscSectionDestroy(&osec));
343907a3e9cSStefano Zampini
344907a3e9cSStefano Zampini /* material parameters */
345e4d5475eSStefano Zampini /* TODO: what if these values changes? add to some DM hook? */
346e4d5475eSStefano Zampini PetscCall(DMTransferMaterialParameters(dm, sf, odm));
347907a3e9cSStefano Zampini PetscCall(PetscSFDestroy(&sf));
348907a3e9cSStefano Zampini
349907a3e9cSStefano Zampini PetscCall(DMViewFromOptions(dm, NULL, "-dm_plex_view_neumann_original"));
350907a3e9cSStefano Zampini PetscCall(PetscObjectSetName((PetscObject)odm, "OVL"));
351907a3e9cSStefano Zampini PetscCall(DMViewFromOptions(odm, NULL, "-dm_plex_view_neumann_overlap"));
352907a3e9cSStefano Zampini
353907a3e9cSStefano Zampini /* MATIS for the overlap region
354907a3e9cSStefano Zampini the HPDDM interface wants local matrices,
355907a3e9cSStefano Zampini we attach the global MATIS to the overlap IS,
356907a3e9cSStefano Zampini since we need it to do assembly */
357907a3e9cSStefano Zampini PetscCall(DMSetMatType(odm, MATIS));
358907a3e9cSStefano Zampini PetscCall(DMCreateMatrix(odm, &pJ));
359907a3e9cSStefano Zampini PetscCall(MatISGetLocalMat(pJ, J));
360907a3e9cSStefano Zampini PetscCall(PetscObjectReference((PetscObject)*J));
361907a3e9cSStefano Zampini
362907a3e9cSStefano Zampini /* overlap IS */
363907a3e9cSStefano Zampini PetscCall(MatISGetLocalToGlobalMapping(pJ, &l2g, NULL));
364907a3e9cSStefano Zampini PetscCall(ISLocalToGlobalMappingGetSize(l2g, &n));
365907a3e9cSStefano Zampini PetscCall(ISLocalToGlobalMappingGetIndices(l2g, &idxs));
366907a3e9cSStefano Zampini PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), n, idxs, PETSC_COPY_VALUES, ovl));
367907a3e9cSStefano Zampini PetscCall(ISLocalToGlobalMappingRestoreIndices(l2g, &idxs));
368907a3e9cSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)*ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject)pJ));
369907a3e9cSStefano Zampini PetscCall(DMDestroy(&odm));
370907a3e9cSStefano Zampini PetscCall(MatDestroy(&pJ));
371907a3e9cSStefano Zampini
372907a3e9cSStefano Zampini /* special purpose setup function (composed in DMPlexSetSNESLocalFEM) */
373907a3e9cSStefano Zampini PetscCall(PetscObjectQueryFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", setup));
374907a3e9cSStefano Zampini if (*setup) PetscCall(PetscObjectCompose((PetscObject)*ovl, "_DM_Original_HPDDM", (PetscObject)dm));
375907a3e9cSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
376907a3e9cSStefano Zampini }
377