xref: /petsc/src/dm/impls/plex/plexdd.c (revision e4d5475ea39a6400e1a24d6c1cc2000cd57a112d)
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 
6*e4d5475eSStefano Zampini static PetscErrorCode DMTransferMaterialParameters(DM dm, PetscSF sf, DM odm)
7*e4d5475eSStefano Zampini {
8*e4d5475eSStefano Zampini   Vec A;
9*e4d5475eSStefano Zampini 
10*e4d5475eSStefano Zampini   PetscFunctionBegin;
11*e4d5475eSStefano Zampini   /* TODO handle regions? */
12*e4d5475eSStefano Zampini   PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &A));
13*e4d5475eSStefano Zampini   if (A) {
14*e4d5475eSStefano Zampini     DM           dmAux, ocdm, odmAux;
15*e4d5475eSStefano Zampini     Vec          oA, oAt;
16*e4d5475eSStefano Zampini     PetscSection sec, osec;
17*e4d5475eSStefano Zampini 
18*e4d5475eSStefano Zampini     PetscCall(VecGetDM(A, &dmAux));
19*e4d5475eSStefano Zampini     PetscCall(DMClone(odm, &odmAux));
20*e4d5475eSStefano Zampini     PetscCall(DMGetCoordinateDM(odm, &ocdm));
21*e4d5475eSStefano Zampini     PetscCall(DMSetCoordinateDM(odmAux, ocdm));
22*e4d5475eSStefano Zampini     PetscCall(DMCopyDisc(dmAux, odmAux));
23*e4d5475eSStefano Zampini 
24*e4d5475eSStefano Zampini     PetscCall(DMGetLocalSection(dmAux, &sec));
25*e4d5475eSStefano Zampini     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)sec), &osec));
26*e4d5475eSStefano Zampini     PetscCall(VecCreate(PetscObjectComm((PetscObject)A), &oAt));
27*e4d5475eSStefano Zampini     PetscCall(DMPlexDistributeField(dmAux, sf, sec, A, osec, oAt));
28*e4d5475eSStefano Zampini     PetscCall(DMSetLocalSection(odmAux, osec));
29*e4d5475eSStefano Zampini     PetscCall(PetscSectionDestroy(&osec));
30*e4d5475eSStefano Zampini     PetscCall(DMCreateLocalVector(odmAux, &oA));
31*e4d5475eSStefano Zampini     PetscCall(DMDestroy(&odmAux));
32*e4d5475eSStefano Zampini     PetscCall(VecCopy(oAt, oA));
33*e4d5475eSStefano Zampini     PetscCall(VecDestroy(&oAt));
34*e4d5475eSStefano Zampini     PetscCall(DMSetAuxiliaryVec(odm, NULL, 0, 0, oA));
35*e4d5475eSStefano Zampini     PetscCall(VecDestroy(&oA));
36*e4d5475eSStefano Zampini   }
37*e4d5475eSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
38*e4d5475eSStefano Zampini }
39*e4d5475eSStefano 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   }
91*e4d5475eSStefano Zampini   /* TODO: what if these values changes? add to some DM hook? */
92*e4d5475eSStefano Zampini   PetscCall(DMTransferMaterialParameters(dm, migrationSF, odm));
93907a3e9cSStefano Zampini 
94907a3e9cSStefano Zampini   /* TODO: it would be nice to automatically translate filenames with wildcards:
95907a3e9cSStefano Zampini            name%r.vtu -> name${rank}.vtu
96907a3e9cSStefano Zampini            name%R.vtu -> name${PetscGlobalRank}.vtu
97907a3e9cSStefano Zampini            So that -dm_view vtk:name%R.vtu will automatically translate to the commented code below
98907a3e9cSStefano Zampini   */
99907a3e9cSStefano Zampini #if 0
100907a3e9cSStefano Zampini   {
101907a3e9cSStefano Zampini     char        name[256];
102907a3e9cSStefano Zampini     PetscViewer viewer;
103907a3e9cSStefano Zampini 
104907a3e9cSStefano Zampini     PetscCall(PetscSNPrintf(name, sizeof(name), "ovl_mesh_%d.vtu", PetscGlobalRank));
105907a3e9cSStefano Zampini     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)odm), name, FILE_MODE_WRITE, &viewer));
106907a3e9cSStefano Zampini     PetscCall(DMView(odm, viewer));
107907a3e9cSStefano Zampini     PetscCall(PetscViewerDestroy(&viewer));
108907a3e9cSStefano Zampini   }
109907a3e9cSStefano Zampini #endif
110907a3e9cSStefano Zampini   PetscCall(DMViewFromOptions(odm, (PetscObject)dm, "-dm_plex_dd_overlap_dm_view"));
111907a3e9cSStefano Zampini 
112907a3e9cSStefano Zampini   /* propagate original global ordering to overlapping DM */
113907a3e9cSStefano Zampini   PetscCall(DMGetSectionSF(dm, &sectionSF));
114907a3e9cSStefano Zampini   PetscCall(DMGetLocalSection(dm, &sec));
115907a3e9cSStefano Zampini   PetscCall(PetscSectionGetStorageSize(sec, &nl));
116907a3e9cSStefano Zampini   PetscCall(DMGetGlobalVector(dm, &gvec));
117907a3e9cSStefano Zampini   PetscCall(VecGetOwnershipRange(gvec, &rst, &ren));
118907a3e9cSStefano Zampini   PetscCall(DMRestoreGlobalVector(dm, &gvec));
119907a3e9cSStefano Zampini   PetscCall(ISCreateStride(PETSC_COMM_SELF, ren - rst, rst, 1, &gi_is)); /* non-overlapping dofs */
120907a3e9cSStefano Zampini   PetscCall(PetscMalloc1(nl, &lidxs));
121907a3e9cSStefano Zampini   for (PetscInt i = 0; i < nl; i++) lidxs[i] = -1;
122907a3e9cSStefano Zampini   PetscCall(ISGetIndices(gi_is, (const PetscInt **)&gidxs));
123907a3e9cSStefano Zampini   PetscCall(PetscSFBcastBegin(sectionSF, MPIU_INT, gidxs, lidxs, MPI_REPLACE));
124907a3e9cSStefano Zampini   PetscCall(PetscSFBcastEnd(sectionSF, MPIU_INT, gidxs, lidxs, MPI_REPLACE));
125907a3e9cSStefano Zampini   PetscCall(ISRestoreIndices(gi_is, (const PetscInt **)&gidxs));
126907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), nl, lidxs, PETSC_OWN_POINTER, &lis));
127907a3e9cSStefano Zampini   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)odm), &tsec));
128907a3e9cSStefano Zampini   PetscCall(DMPlexDistributeFieldIS(dm, migrationSF, sec, lis, tsec, &gis));
129907a3e9cSStefano Zampini   PetscCall(PetscSectionDestroy(&tsec));
130907a3e9cSStefano Zampini   PetscCall(ISDestroy(&lis));
131907a3e9cSStefano Zampini   PetscCall(PetscSFDestroy(&migrationSF));
132907a3e9cSStefano Zampini 
133907a3e9cSStefano Zampini   /* make dofs on the overlap boundary (not the global boundary) constrained */
134907a3e9cSStefano Zampini   PetscCall(DMGetLabel(odm, oname, &label));
135907a3e9cSStefano Zampini   PetscCall(DMPlexLabelComplete(odm, label));
136907a3e9cSStefano Zampini   PetscCall(DMGetLocalSection(odm, &tsec));
137907a3e9cSStefano Zampini   PetscCall(PetscSectionGetChart(tsec, &pStart, &pEnd));
138907a3e9cSStefano Zampini   PetscCall(DMLabelCreateIndex(label, pStart, pEnd));
139907a3e9cSStefano Zampini   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)tsec), &sec));
140907a3e9cSStefano Zampini   PetscCall(PetscSectionCopy_Internal(tsec, sec, label->bt));
141907a3e9cSStefano Zampini   PetscCall(DMSetLocalSection(odm, sec));
142907a3e9cSStefano Zampini   PetscCall(PetscSectionDestroy(&sec));
143907a3e9cSStefano Zampini   PetscCall(DMRemoveLabel(odm, oname, NULL));
144907a3e9cSStefano Zampini 
145907a3e9cSStefano Zampini   /* Create index sets for dofs in the overlap dm */
146907a3e9cSStefano Zampini   PetscCall(DMGetSectionSF(odm, &sectionSF));
147907a3e9cSStefano Zampini   PetscCall(DMGetLocalSection(odm, &olsec));
148907a3e9cSStefano Zampini   PetscCall(DMGetGlobalSection(odm, &ogsec));
149907a3e9cSStefano Zampini   PetscCall(PetscSectionViewFromOptions(ogsec, (PetscObject)dm, "-dm_plex_dd_overlap_gsection_view"));
150907a3e9cSStefano Zampini   PetscCall(PetscSectionViewFromOptions(olsec, (PetscObject)dm, "-dm_plex_dd_overlap_lsection_view"));
151907a3e9cSStefano Zampini   ni = ren - rst;
152907a3e9cSStefano Zampini   PetscCall(PetscSectionGetConstrainedStorageSize(ogsec, &no)); /* dofs in the overlap */
153907a3e9cSStefano Zampini   PetscCall(PetscSectionGetStorageSize(olsec, &nl));            /* local dofs in the overlap */
154907a3e9cSStefano Zampini   PetscCall(PetscMalloc1(no, &gidxs));
155907a3e9cSStefano Zampini   PetscCall(ISGetIndices(gis, (const PetscInt **)&lidxs));
156907a3e9cSStefano Zampini   PetscCall(PetscSFReduceBegin(sectionSF, MPIU_INT, lidxs, gidxs, MPI_REPLACE));
157907a3e9cSStefano Zampini   PetscCall(PetscSFReduceEnd(sectionSF, MPIU_INT, lidxs, gidxs, MPI_REPLACE));
158907a3e9cSStefano Zampini   PetscCall(ISRestoreIndices(gis, (const PetscInt **)&lidxs));
159907a3e9cSStefano Zampini 
160907a3e9cSStefano Zampini   /* non-overlapping dofs */
161907a3e9cSStefano Zampini   PetscCall(PetscMalloc1(no, &lidxs));
162907a3e9cSStefano Zampini   c = 0;
163907a3e9cSStefano Zampini   for (PetscInt i = 0; i < no; i++) {
164907a3e9cSStefano Zampini     if (gidxs[i] >= rst && gidxs[i] < ren) lidxs[c++] = i;
165907a3e9cSStefano Zampini   }
166907a3e9cSStefano Zampini   PetscCheck(c == ni, PETSC_COMM_SELF, PETSC_ERR_PLIB, "%" PetscInt_FMT " != %" PetscInt_FMT "\n", c, ni);
167907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), ni, lidxs, PETSC_OWN_POINTER, &li_is));
168907a3e9cSStefano Zampini 
169907a3e9cSStefano Zampini   /* global dofs in the overlap */
170907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), no, gidxs, PETSC_OWN_POINTER, &go_is));
171907a3e9cSStefano Zampini   PetscCall(ISViewFromOptions(go_is, (PetscObject)dm, "-dm_plex_dd_overlap_gois_view"));
172907a3e9cSStefano Zampini   /* PetscCall(ISCreateStride(PetscObjectComm((PetscObject)odm), no, 0, 1, &lo_is)); */
173907a3e9cSStefano Zampini 
174907a3e9cSStefano Zampini   /* local dofs of the overlapping subdomain (we actually need only dofs on the boundary of the subdomain) */
175907a3e9cSStefano Zampini   PetscCall(PetscMalloc1(nl, &lidxs));
176907a3e9cSStefano Zampini   PetscCall(PetscMalloc1(nl, &gidxs));
177907a3e9cSStefano Zampini   PetscCall(ISGetIndices(gis, (const PetscInt **)&tidxs));
178907a3e9cSStefano Zampini   c = 0;
179907a3e9cSStefano Zampini   for (PetscInt i = 0; i < nl; i++) {
180907a3e9cSStefano Zampini     if (tidxs[i] < 0) continue;
181907a3e9cSStefano Zampini     lidxs[c] = i;
182907a3e9cSStefano Zampini     gidxs[c] = tidxs[i];
183907a3e9cSStefano Zampini     c++;
184907a3e9cSStefano Zampini   }
185907a3e9cSStefano Zampini   PetscCall(ISRestoreIndices(gis, (const PetscInt **)&tidxs));
186907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), c, gidxs, PETSC_OWN_POINTER, &gl_is));
187907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), c, lidxs, PETSC_OWN_POINTER, &ll_is));
188907a3e9cSStefano Zampini   PetscCall(ISViewFromOptions(gl_is, (PetscObject)dm, "-dm_plex_dd_overlap_glis_view"));
189907a3e9cSStefano Zampini 
190907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_gi", (PetscObject)gi_is));
191907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_li", (PetscObject)li_is));
192907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_go", (PetscObject)go_is));
193907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_gl", (PetscObject)gl_is));
194907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_ll", (PetscObject)ll_is));
195907a3e9cSStefano Zampini 
196907a3e9cSStefano Zampini   if (innerises) (*innerises)[0] = gi_is;
197907a3e9cSStefano Zampini   else PetscCall(ISDestroy(&gi_is));
198907a3e9cSStefano Zampini   if (outerises) (*outerises)[0] = go_is;
199907a3e9cSStefano Zampini   else PetscCall(ISDestroy(&go_is));
200907a3e9cSStefano Zampini   if (dms) (*dms)[0] = odm;
201907a3e9cSStefano Zampini   else PetscCall(DMDestroy(&odm));
202907a3e9cSStefano Zampini   PetscCall(ISDestroy(&li_is));
203907a3e9cSStefano Zampini   PetscCall(ISDestroy(&gl_is));
204907a3e9cSStefano Zampini   PetscCall(ISDestroy(&ll_is));
205907a3e9cSStefano Zampini   PetscCall(ISDestroy(&gis));
206907a3e9cSStefano Zampini 
207907a3e9cSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
208907a3e9cSStefano Zampini }
209907a3e9cSStefano Zampini 
210907a3e9cSStefano Zampini PetscErrorCode DMCreateDomainDecompositionScatters_Plex(DM dm, PetscInt n, DM *subdms, VecScatter **iscat, VecScatter **oscat, VecScatter **lscat)
211907a3e9cSStefano Zampini {
212907a3e9cSStefano Zampini   Vec gvec, svec, lvec;
213*e4d5475eSStefano Zampini   IS  gi_is, li_is, go_is, gl_is, ll_is;
214907a3e9cSStefano Zampini 
215907a3e9cSStefano Zampini   PetscFunctionBegin;
216907a3e9cSStefano Zampini   if (iscat) PetscCall(PetscMalloc1(n, iscat));
217907a3e9cSStefano Zampini   if (oscat) PetscCall(PetscMalloc1(n, oscat));
218907a3e9cSStefano Zampini   if (lscat) PetscCall(PetscMalloc1(n, lscat));
219907a3e9cSStefano Zampini 
220907a3e9cSStefano Zampini   PetscCall(DMGetGlobalVector(dm, &gvec));
221907a3e9cSStefano Zampini   for (PetscInt i = 0; i < n; i++) {
222907a3e9cSStefano Zampini     PetscCall(DMGetGlobalVector(subdms[i], &svec));
223907a3e9cSStefano Zampini     PetscCall(DMGetLocalVector(subdms[i], &lvec));
224907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_gi", (PetscObject *)&gi_is));
225907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_li", (PetscObject *)&li_is));
226907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_go", (PetscObject *)&go_is));
227907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_gl", (PetscObject *)&gl_is));
228907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_ll", (PetscObject *)&ll_is));
229*e4d5475eSStefano Zampini     PetscCheck(gi_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition");
230*e4d5475eSStefano Zampini     PetscCheck(li_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition");
231*e4d5475eSStefano Zampini     PetscCheck(go_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition");
232*e4d5475eSStefano Zampini     PetscCheck(gl_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition");
233*e4d5475eSStefano Zampini     PetscCheck(ll_is, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "SubDM not obtained from DMCreateDomainDecomposition");
234907a3e9cSStefano Zampini     if (iscat) PetscCall(VecScatterCreate(gvec, gi_is, svec, li_is, &(*iscat)[i]));
235*e4d5475eSStefano Zampini     if (oscat) PetscCall(VecScatterCreate(gvec, go_is, svec, NULL, &(*oscat)[i]));
236907a3e9cSStefano Zampini     if (lscat) PetscCall(VecScatterCreate(gvec, gl_is, lvec, ll_is, &(*lscat)[i]));
237907a3e9cSStefano Zampini     PetscCall(DMRestoreGlobalVector(subdms[i], &svec));
238907a3e9cSStefano Zampini     PetscCall(DMRestoreLocalVector(subdms[i], &lvec));
239907a3e9cSStefano Zampini   }
240907a3e9cSStefano Zampini   PetscCall(DMRestoreGlobalVector(dm, &gvec));
241907a3e9cSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
242907a3e9cSStefano Zampini }
243907a3e9cSStefano Zampini 
244907a3e9cSStefano Zampini /*
245907a3e9cSStefano Zampini      DMCreateNeumannOverlap_Plex - Generates an IS, an unassembled (Neumann) Mat, a setup function, and the corresponding context to be used by PCHPDDM.
246907a3e9cSStefano Zampini 
247907a3e9cSStefano Zampini    Input Parameter:
248907a3e9cSStefano Zampini .     dm - preconditioner context
249907a3e9cSStefano Zampini 
250907a3e9cSStefano Zampini    Output Parameters:
251907a3e9cSStefano Zampini +     ovl - index set of overlapping subdomains
252907a3e9cSStefano Zampini .     J - unassembled (Neumann) local matrix
253907a3e9cSStefano Zampini .     setup - function for generating the matrix
254907a3e9cSStefano Zampini -     setup_ctx - context for setup
255907a3e9cSStefano Zampini 
256907a3e9cSStefano Zampini    Options Database Keys:
257907a3e9cSStefano Zampini +   -dm_plex_view_neumann_original - view the DM without overlap
258907a3e9cSStefano Zampini -   -dm_plex_view_neumann_overlap - view the DM with overlap as needed by PCHPDDM
259907a3e9cSStefano Zampini 
260907a3e9cSStefano Zampini    Level: advanced
261907a3e9cSStefano Zampini 
262907a3e9cSStefano Zampini .seealso: `DMCreate()`, `DM`, `MATIS`, `PCHPDDM`, `PCHPDDMSetAuxiliaryMat()`
263907a3e9cSStefano Zampini */
264907a3e9cSStefano Zampini PetscErrorCode DMCreateNeumannOverlap_Plex(DM dm, IS *ovl, Mat *J, PetscErrorCode (**setup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void **setup_ctx)
265907a3e9cSStefano Zampini {
266907a3e9cSStefano Zampini   DM                     odm;
267907a3e9cSStefano Zampini   Mat                    pJ;
268907a3e9cSStefano Zampini   PetscSF                sf = NULL;
269907a3e9cSStefano Zampini   PetscSection           sec, osec;
270907a3e9cSStefano Zampini   ISLocalToGlobalMapping l2g;
271907a3e9cSStefano Zampini   const PetscInt        *idxs;
272907a3e9cSStefano Zampini   PetscInt               n, mh;
273907a3e9cSStefano Zampini 
274907a3e9cSStefano Zampini   PetscFunctionBegin;
275907a3e9cSStefano Zampini   *setup     = NULL;
276907a3e9cSStefano Zampini   *setup_ctx = NULL;
277907a3e9cSStefano Zampini   *ovl       = NULL;
278907a3e9cSStefano Zampini   *J         = NULL;
279907a3e9cSStefano Zampini 
280907a3e9cSStefano Zampini   /* Overlapped mesh
281907a3e9cSStefano Zampini      overlap is a little more generous, since it is not computed starting from the owned (Dirichlet) points, but from the locally owned cells */
282907a3e9cSStefano Zampini   PetscCall(DMPlexDistributeOverlap(dm, 1, &sf, &odm));
283907a3e9cSStefano Zampini   if (!odm) {
284907a3e9cSStefano Zampini     PetscCall(PetscSFDestroy(&sf));
285907a3e9cSStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
286907a3e9cSStefano Zampini   }
287907a3e9cSStefano Zampini 
288907a3e9cSStefano Zampini   /* share discretization */
289907a3e9cSStefano Zampini   PetscCall(DMGetLocalSection(dm, &sec));
290907a3e9cSStefano Zampini   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)sec), &osec));
291907a3e9cSStefano Zampini   PetscCall(PetscSFDistributeSection(sf, sec, NULL, osec));
292907a3e9cSStefano Zampini   /* what to do here? using both is fine? */
293907a3e9cSStefano Zampini   PetscCall(DMSetLocalSection(odm, osec));
294907a3e9cSStefano Zampini   PetscCall(DMCopyDisc(dm, odm));
295907a3e9cSStefano Zampini   PetscCall(DMPlexGetMaxProjectionHeight(dm, &mh));
296907a3e9cSStefano Zampini   PetscCall(DMPlexSetMaxProjectionHeight(odm, mh));
297907a3e9cSStefano Zampini   PetscCall(PetscSectionDestroy(&osec));
298907a3e9cSStefano Zampini 
299907a3e9cSStefano Zampini   /* material parameters */
300*e4d5475eSStefano Zampini   /* TODO: what if these values changes? add to some DM hook? */
301*e4d5475eSStefano Zampini   PetscCall(DMTransferMaterialParameters(dm, sf, odm));
302907a3e9cSStefano Zampini   PetscCall(PetscSFDestroy(&sf));
303907a3e9cSStefano Zampini 
304907a3e9cSStefano Zampini   PetscCall(DMViewFromOptions(dm, NULL, "-dm_plex_view_neumann_original"));
305907a3e9cSStefano Zampini   PetscCall(PetscObjectSetName((PetscObject)odm, "OVL"));
306907a3e9cSStefano Zampini   PetscCall(DMViewFromOptions(odm, NULL, "-dm_plex_view_neumann_overlap"));
307907a3e9cSStefano Zampini 
308907a3e9cSStefano Zampini   /* MATIS for the overlap region
309907a3e9cSStefano Zampini      the HPDDM interface wants local matrices,
310907a3e9cSStefano Zampini      we attach the global MATIS to the overlap IS,
311907a3e9cSStefano Zampini      since we need it to do assembly */
312907a3e9cSStefano Zampini   PetscCall(DMSetMatType(odm, MATIS));
313907a3e9cSStefano Zampini   PetscCall(DMCreateMatrix(odm, &pJ));
314907a3e9cSStefano Zampini   PetscCall(MatISGetLocalMat(pJ, J));
315907a3e9cSStefano Zampini   PetscCall(PetscObjectReference((PetscObject)*J));
316907a3e9cSStefano Zampini 
317907a3e9cSStefano Zampini   /* overlap IS */
318907a3e9cSStefano Zampini   PetscCall(MatISGetLocalToGlobalMapping(pJ, &l2g, NULL));
319907a3e9cSStefano Zampini   PetscCall(ISLocalToGlobalMappingGetSize(l2g, &n));
320907a3e9cSStefano Zampini   PetscCall(ISLocalToGlobalMappingGetIndices(l2g, &idxs));
321907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), n, idxs, PETSC_COPY_VALUES, ovl));
322907a3e9cSStefano Zampini   PetscCall(ISLocalToGlobalMappingRestoreIndices(l2g, &idxs));
323907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)*ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject)pJ));
324907a3e9cSStefano Zampini   PetscCall(DMDestroy(&odm));
325907a3e9cSStefano Zampini   PetscCall(MatDestroy(&pJ));
326907a3e9cSStefano Zampini 
327907a3e9cSStefano Zampini   /* special purpose setup function (composed in DMPlexSetSNESLocalFEM) */
328907a3e9cSStefano Zampini   PetscCall(PetscObjectQueryFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", setup));
329907a3e9cSStefano Zampini   if (*setup) PetscCall(PetscObjectCompose((PetscObject)*ovl, "_DM_Original_HPDDM", (PetscObject)dm));
330907a3e9cSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
331907a3e9cSStefano Zampini }
332