xref: /petsc/src/dm/impls/plex/plexdd.c (revision 907a3e9c7e62a1244ef44b3544cf5c6924719401)
1*907a3e9cSStefano Zampini #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
2*907a3e9cSStefano Zampini #include <petsc/private/dmlabelimpl.h>
3*907a3e9cSStefano Zampini #include <petsc/private/sectionimpl.h>
4*907a3e9cSStefano Zampini #include <petsc/private/sfimpl.h>
5*907a3e9cSStefano Zampini 
6*907a3e9cSStefano Zampini PetscErrorCode DMCreateDomainDecomposition_Plex(DM dm, PetscInt *nsub, char ***names, IS **innerises, IS **outerises, DM **dms)
7*907a3e9cSStefano Zampini {
8*907a3e9cSStefano Zampini   DM           odm;
9*907a3e9cSStefano Zampini   PetscSF      migrationSF = NULL, sectionSF;
10*907a3e9cSStefano Zampini   PetscSection sec, tsec, ogsec, olsec;
11*907a3e9cSStefano Zampini   PetscInt     n, mh, ddovl = 0, pStart, pEnd, ni, no, nl;
12*907a3e9cSStefano Zampini   PetscDS      ds;
13*907a3e9cSStefano Zampini   DMLabel      label;
14*907a3e9cSStefano Zampini   const char  *oname = "__internal_plex_dd_ovl_";
15*907a3e9cSStefano Zampini   IS           gi_is, li_is, go_is, gl_is, ll_is;
16*907a3e9cSStefano Zampini   IS           gis, lis;
17*907a3e9cSStefano Zampini   PetscInt     rst, ren, c, *gidxs, *lidxs, *tidxs;
18*907a3e9cSStefano Zampini   Vec          gvec;
19*907a3e9cSStefano Zampini 
20*907a3e9cSStefano Zampini   PetscFunctionBegin;
21*907a3e9cSStefano Zampini   n = 1;
22*907a3e9cSStefano Zampini   if (nsub) *nsub = n;
23*907a3e9cSStefano Zampini   if (names) PetscCall(PetscCalloc1(n, names));
24*907a3e9cSStefano Zampini   if (innerises) PetscCall(PetscCalloc1(n, innerises));
25*907a3e9cSStefano Zampini   if (outerises) PetscCall(PetscCalloc1(n, outerises));
26*907a3e9cSStefano Zampini   if (dms) PetscCall(PetscCalloc1(n, dms));
27*907a3e9cSStefano Zampini 
28*907a3e9cSStefano Zampini   PetscObjectOptionsBegin((PetscObject)dm);
29*907a3e9cSStefano Zampini   PetscCall(PetscOptionsBoundedInt("-dm_plex_dd_overlap", "The size of the overlap halo for the subdomains", "DMCreateDomainDecomposition", ddovl, &ddovl, NULL, 0));
30*907a3e9cSStefano Zampini   PetscOptionsEnd();
31*907a3e9cSStefano Zampini 
32*907a3e9cSStefano Zampini   PetscCall(DMViewFromOptions(dm, NULL, "-dm_plex_dd_dm_view"));
33*907a3e9cSStefano Zampini   PetscCall(DMPlexDistributeOverlap_Internal(dm, ddovl + 1, PETSC_COMM_SELF, oname, &migrationSF, &odm));
34*907a3e9cSStefano Zampini   if (!odm) PetscCall(DMClone(dm, &odm));
35*907a3e9cSStefano Zampini   if (migrationSF) PetscCall(PetscSFViewFromOptions(migrationSF, (PetscObject)dm, "-dm_plex_dd_sf_view"));
36*907a3e9cSStefano Zampini 
37*907a3e9cSStefano Zampini   PetscCall(DMPlexGetMaxProjectionHeight(dm, &mh));
38*907a3e9cSStefano Zampini   PetscCall(DMPlexSetMaxProjectionHeight(odm, mh));
39*907a3e9cSStefano Zampini 
40*907a3e9cSStefano Zampini   /* share discretization */
41*907a3e9cSStefano Zampini   /* TODO material parameters */
42*907a3e9cSStefano Zampini   /* TODO Labels for regions may need to updated,
43*907a3e9cSStefano Zampini      now it uses the original ones, not the ones for the odm.
44*907a3e9cSStefano Zampini      Not sure what to do here */
45*907a3e9cSStefano Zampini   /* PetscCall(DMCopyDisc(dm, odm)); */
46*907a3e9cSStefano Zampini   PetscCall(DMGetDS(odm, &ds));
47*907a3e9cSStefano Zampini   if (!ds) {
48*907a3e9cSStefano Zampini     PetscCall(DMGetLocalSection(dm, &sec));
49*907a3e9cSStefano Zampini     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)odm), &tsec));
50*907a3e9cSStefano Zampini     if (migrationSF) {
51*907a3e9cSStefano Zampini       PetscCall(PetscSFDistributeSection(migrationSF, sec, NULL, tsec));
52*907a3e9cSStefano Zampini     } else {
53*907a3e9cSStefano Zampini       PetscCall(PetscSectionCopy(sec, tsec));
54*907a3e9cSStefano Zampini     }
55*907a3e9cSStefano Zampini     PetscCall(DMSetLocalSection(dm, tsec));
56*907a3e9cSStefano Zampini     PetscCall(PetscSectionDestroy(&tsec));
57*907a3e9cSStefano Zampini   }
58*907a3e9cSStefano Zampini 
59*907a3e9cSStefano Zampini   /* TODO: it would be nice to automatically translate filenames with wildcards:
60*907a3e9cSStefano Zampini            name%r.vtu -> name${rank}.vtu
61*907a3e9cSStefano Zampini            name%R.vtu -> name${PetscGlobalRank}.vtu
62*907a3e9cSStefano Zampini            So that -dm_view vtk:name%R.vtu will automatically translate to the commented code below
63*907a3e9cSStefano Zampini   */
64*907a3e9cSStefano Zampini #if 0
65*907a3e9cSStefano Zampini   {
66*907a3e9cSStefano Zampini     char        name[256];
67*907a3e9cSStefano Zampini     PetscViewer viewer;
68*907a3e9cSStefano Zampini 
69*907a3e9cSStefano Zampini     PetscCall(PetscSNPrintf(name, sizeof(name), "ovl_mesh_%d.vtu", PetscGlobalRank));
70*907a3e9cSStefano Zampini     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)odm), name, FILE_MODE_WRITE, &viewer));
71*907a3e9cSStefano Zampini     PetscCall(DMView(odm, viewer));
72*907a3e9cSStefano Zampini     PetscCall(PetscViewerDestroy(&viewer));
73*907a3e9cSStefano Zampini   }
74*907a3e9cSStefano Zampini #endif
75*907a3e9cSStefano Zampini   PetscCall(DMViewFromOptions(odm, (PetscObject)dm, "-dm_plex_dd_overlap_dm_view"));
76*907a3e9cSStefano Zampini 
77*907a3e9cSStefano Zampini   /* propagate original global ordering to overlapping DM */
78*907a3e9cSStefano Zampini   PetscCall(DMGetSectionSF(dm, &sectionSF));
79*907a3e9cSStefano Zampini   PetscCall(DMGetLocalSection(dm, &sec));
80*907a3e9cSStefano Zampini   PetscCall(PetscSectionGetStorageSize(sec, &nl));
81*907a3e9cSStefano Zampini   PetscCall(DMGetGlobalVector(dm, &gvec));
82*907a3e9cSStefano Zampini   PetscCall(VecGetOwnershipRange(gvec, &rst, &ren));
83*907a3e9cSStefano Zampini   PetscCall(DMRestoreGlobalVector(dm, &gvec));
84*907a3e9cSStefano Zampini   PetscCall(ISCreateStride(PETSC_COMM_SELF, ren - rst, rst, 1, &gi_is)); /* non-overlapping dofs */
85*907a3e9cSStefano Zampini   PetscCall(PetscMalloc1(nl, &lidxs));
86*907a3e9cSStefano Zampini   for (PetscInt i = 0; i < nl; i++) lidxs[i] = -1;
87*907a3e9cSStefano Zampini   PetscCall(ISGetIndices(gi_is, (const PetscInt **)&gidxs));
88*907a3e9cSStefano Zampini   PetscCall(PetscSFBcastBegin(sectionSF, MPIU_INT, gidxs, lidxs, MPI_REPLACE));
89*907a3e9cSStefano Zampini   PetscCall(PetscSFBcastEnd(sectionSF, MPIU_INT, gidxs, lidxs, MPI_REPLACE));
90*907a3e9cSStefano Zampini   PetscCall(ISRestoreIndices(gi_is, (const PetscInt **)&gidxs));
91*907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), nl, lidxs, PETSC_OWN_POINTER, &lis));
92*907a3e9cSStefano Zampini   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)odm), &tsec));
93*907a3e9cSStefano Zampini   PetscCall(DMPlexDistributeFieldIS(dm, migrationSF, sec, lis, tsec, &gis));
94*907a3e9cSStefano Zampini   PetscCall(PetscSectionDestroy(&tsec));
95*907a3e9cSStefano Zampini   PetscCall(ISDestroy(&lis));
96*907a3e9cSStefano Zampini   PetscCall(PetscSFDestroy(&migrationSF));
97*907a3e9cSStefano Zampini 
98*907a3e9cSStefano Zampini   /* make dofs on the overlap boundary (not the global boundary) constrained */
99*907a3e9cSStefano Zampini   PetscCall(DMGetLabel(odm, oname, &label));
100*907a3e9cSStefano Zampini   PetscCall(DMPlexLabelComplete(odm, label));
101*907a3e9cSStefano Zampini   PetscCall(DMGetLocalSection(odm, &tsec));
102*907a3e9cSStefano Zampini   PetscCall(PetscSectionGetChart(tsec, &pStart, &pEnd));
103*907a3e9cSStefano Zampini   PetscCall(DMLabelCreateIndex(label, pStart, pEnd));
104*907a3e9cSStefano Zampini   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)tsec), &sec));
105*907a3e9cSStefano Zampini   PetscCall(PetscSectionCopy_Internal(tsec, sec, label->bt));
106*907a3e9cSStefano Zampini   PetscCall(DMSetLocalSection(odm, sec));
107*907a3e9cSStefano Zampini   PetscCall(PetscSectionDestroy(&sec));
108*907a3e9cSStefano Zampini   PetscCall(DMRemoveLabel(odm, oname, NULL));
109*907a3e9cSStefano Zampini 
110*907a3e9cSStefano Zampini   /* Create index sets for dofs in the overlap dm */
111*907a3e9cSStefano Zampini   PetscCall(DMGetSectionSF(odm, &sectionSF));
112*907a3e9cSStefano Zampini   PetscCall(DMGetLocalSection(odm, &olsec));
113*907a3e9cSStefano Zampini   PetscCall(DMGetGlobalSection(odm, &ogsec));
114*907a3e9cSStefano Zampini   PetscCall(PetscSectionViewFromOptions(ogsec, (PetscObject)dm, "-dm_plex_dd_overlap_gsection_view"));
115*907a3e9cSStefano Zampini   PetscCall(PetscSectionViewFromOptions(olsec, (PetscObject)dm, "-dm_plex_dd_overlap_lsection_view"));
116*907a3e9cSStefano Zampini   ni = ren - rst;
117*907a3e9cSStefano Zampini   PetscCall(PetscSectionGetConstrainedStorageSize(ogsec, &no)); /* dofs in the overlap */
118*907a3e9cSStefano Zampini   PetscCall(PetscSectionGetStorageSize(olsec, &nl));            /* local dofs in the overlap */
119*907a3e9cSStefano Zampini   PetscCall(PetscMalloc1(no, &gidxs));
120*907a3e9cSStefano Zampini   PetscCall(ISGetIndices(gis, (const PetscInt **)&lidxs));
121*907a3e9cSStefano Zampini   PetscCall(PetscSFReduceBegin(sectionSF, MPIU_INT, lidxs, gidxs, MPI_REPLACE));
122*907a3e9cSStefano Zampini   PetscCall(PetscSFReduceEnd(sectionSF, MPIU_INT, lidxs, gidxs, MPI_REPLACE));
123*907a3e9cSStefano Zampini   PetscCall(ISRestoreIndices(gis, (const PetscInt **)&lidxs));
124*907a3e9cSStefano Zampini 
125*907a3e9cSStefano Zampini   /* non-overlapping dofs */
126*907a3e9cSStefano Zampini   PetscCall(PetscMalloc1(no, &lidxs));
127*907a3e9cSStefano Zampini   c = 0;
128*907a3e9cSStefano Zampini   for (PetscInt i = 0; i < no; i++) {
129*907a3e9cSStefano Zampini     if (gidxs[i] >= rst && gidxs[i] < ren) lidxs[c++] = i;
130*907a3e9cSStefano Zampini   }
131*907a3e9cSStefano Zampini   PetscCheck(c == ni, PETSC_COMM_SELF, PETSC_ERR_PLIB, "%" PetscInt_FMT " != %" PetscInt_FMT "\n", c, ni);
132*907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), ni, lidxs, PETSC_OWN_POINTER, &li_is));
133*907a3e9cSStefano Zampini 
134*907a3e9cSStefano Zampini   /* global dofs in the overlap */
135*907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), no, gidxs, PETSC_OWN_POINTER, &go_is));
136*907a3e9cSStefano Zampini   PetscCall(ISViewFromOptions(go_is, (PetscObject)dm, "-dm_plex_dd_overlap_gois_view"));
137*907a3e9cSStefano Zampini   /* PetscCall(ISCreateStride(PetscObjectComm((PetscObject)odm), no, 0, 1, &lo_is)); */
138*907a3e9cSStefano Zampini 
139*907a3e9cSStefano Zampini   /* local dofs of the overlapping subdomain (we actually need only dofs on the boundary of the subdomain) */
140*907a3e9cSStefano Zampini   PetscCall(PetscMalloc1(nl, &lidxs));
141*907a3e9cSStefano Zampini   PetscCall(PetscMalloc1(nl, &gidxs));
142*907a3e9cSStefano Zampini   PetscCall(ISGetIndices(gis, (const PetscInt **)&tidxs));
143*907a3e9cSStefano Zampini   c = 0;
144*907a3e9cSStefano Zampini   for (PetscInt i = 0; i < nl; i++) {
145*907a3e9cSStefano Zampini     if (tidxs[i] < 0) continue;
146*907a3e9cSStefano Zampini     lidxs[c] = i;
147*907a3e9cSStefano Zampini     gidxs[c] = tidxs[i];
148*907a3e9cSStefano Zampini     c++;
149*907a3e9cSStefano Zampini   }
150*907a3e9cSStefano Zampini   PetscCall(ISRestoreIndices(gis, (const PetscInt **)&tidxs));
151*907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), c, gidxs, PETSC_OWN_POINTER, &gl_is));
152*907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), c, lidxs, PETSC_OWN_POINTER, &ll_is));
153*907a3e9cSStefano Zampini   PetscCall(ISViewFromOptions(gl_is, (PetscObject)dm, "-dm_plex_dd_overlap_glis_view"));
154*907a3e9cSStefano Zampini 
155*907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_gi", (PetscObject)gi_is));
156*907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_li", (PetscObject)li_is));
157*907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_go", (PetscObject)go_is));
158*907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_lo", NULL));
159*907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_gl", (PetscObject)gl_is));
160*907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)odm, "__Plex_DD_IS_ll", (PetscObject)ll_is));
161*907a3e9cSStefano Zampini 
162*907a3e9cSStefano Zampini   if (innerises) (*innerises)[0] = gi_is;
163*907a3e9cSStefano Zampini   else PetscCall(ISDestroy(&gi_is));
164*907a3e9cSStefano Zampini   if (outerises) (*outerises)[0] = go_is;
165*907a3e9cSStefano Zampini   else PetscCall(ISDestroy(&go_is));
166*907a3e9cSStefano Zampini   if (dms) (*dms)[0] = odm;
167*907a3e9cSStefano Zampini   else PetscCall(DMDestroy(&odm));
168*907a3e9cSStefano Zampini   PetscCall(ISDestroy(&li_is));
169*907a3e9cSStefano Zampini   PetscCall(ISDestroy(&gl_is));
170*907a3e9cSStefano Zampini   PetscCall(ISDestroy(&ll_is));
171*907a3e9cSStefano Zampini   PetscCall(ISDestroy(&gis));
172*907a3e9cSStefano Zampini 
173*907a3e9cSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
174*907a3e9cSStefano Zampini }
175*907a3e9cSStefano Zampini 
176*907a3e9cSStefano Zampini PetscErrorCode DMCreateDomainDecompositionScatters_Plex(DM dm, PetscInt n, DM *subdms, VecScatter **iscat, VecScatter **oscat, VecScatter **lscat)
177*907a3e9cSStefano Zampini {
178*907a3e9cSStefano Zampini   Vec gvec, svec, lvec;
179*907a3e9cSStefano Zampini   IS  gi_is, li_is, go_is, lo_is, gl_is, ll_is;
180*907a3e9cSStefano Zampini 
181*907a3e9cSStefano Zampini   PetscFunctionBegin;
182*907a3e9cSStefano Zampini   if (iscat) PetscCall(PetscMalloc1(n, iscat));
183*907a3e9cSStefano Zampini   if (oscat) PetscCall(PetscMalloc1(n, oscat));
184*907a3e9cSStefano Zampini   if (lscat) PetscCall(PetscMalloc1(n, lscat));
185*907a3e9cSStefano Zampini 
186*907a3e9cSStefano Zampini   PetscCall(DMGetGlobalVector(dm, &gvec));
187*907a3e9cSStefano Zampini   for (PetscInt i = 0; i < n; i++) {
188*907a3e9cSStefano Zampini     PetscCall(DMGetGlobalVector(subdms[i], &svec));
189*907a3e9cSStefano Zampini     PetscCall(DMGetLocalVector(subdms[i], &lvec));
190*907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_gi", (PetscObject *)&gi_is));
191*907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_li", (PetscObject *)&li_is));
192*907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_go", (PetscObject *)&go_is));
193*907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_lo", (PetscObject *)&lo_is));
194*907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_gl", (PetscObject *)&gl_is));
195*907a3e9cSStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)subdms[i], "__Plex_DD_IS_ll", (PetscObject *)&ll_is));
196*907a3e9cSStefano Zampini     if (iscat) PetscCall(VecScatterCreate(gvec, gi_is, svec, li_is, &(*iscat)[i]));
197*907a3e9cSStefano Zampini     if (oscat) PetscCall(VecScatterCreate(gvec, go_is, svec, lo_is, &(*oscat)[i]));
198*907a3e9cSStefano Zampini     if (lscat) PetscCall(VecScatterCreate(gvec, gl_is, lvec, ll_is, &(*lscat)[i]));
199*907a3e9cSStefano Zampini     PetscCall(DMRestoreGlobalVector(subdms[i], &svec));
200*907a3e9cSStefano Zampini     PetscCall(DMRestoreLocalVector(subdms[i], &lvec));
201*907a3e9cSStefano Zampini   }
202*907a3e9cSStefano Zampini   PetscCall(DMRestoreGlobalVector(dm, &gvec));
203*907a3e9cSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
204*907a3e9cSStefano Zampini }
205*907a3e9cSStefano Zampini 
206*907a3e9cSStefano Zampini /*
207*907a3e9cSStefano Zampini      DMCreateNeumannOverlap_Plex - Generates an IS, an unassembled (Neumann) Mat, a setup function, and the corresponding context to be used by PCHPDDM.
208*907a3e9cSStefano Zampini 
209*907a3e9cSStefano Zampini    Input Parameter:
210*907a3e9cSStefano Zampini .     dm - preconditioner context
211*907a3e9cSStefano Zampini 
212*907a3e9cSStefano Zampini    Output Parameters:
213*907a3e9cSStefano Zampini +     ovl - index set of overlapping subdomains
214*907a3e9cSStefano Zampini .     J - unassembled (Neumann) local matrix
215*907a3e9cSStefano Zampini .     setup - function for generating the matrix
216*907a3e9cSStefano Zampini -     setup_ctx - context for setup
217*907a3e9cSStefano Zampini 
218*907a3e9cSStefano Zampini    Options Database Keys:
219*907a3e9cSStefano Zampini +   -dm_plex_view_neumann_original - view the DM without overlap
220*907a3e9cSStefano Zampini -   -dm_plex_view_neumann_overlap - view the DM with overlap as needed by PCHPDDM
221*907a3e9cSStefano Zampini 
222*907a3e9cSStefano Zampini    Level: advanced
223*907a3e9cSStefano Zampini 
224*907a3e9cSStefano Zampini .seealso: `DMCreate()`, `DM`, `MATIS`, `PCHPDDM`, `PCHPDDMSetAuxiliaryMat()`
225*907a3e9cSStefano Zampini */
226*907a3e9cSStefano Zampini PetscErrorCode DMCreateNeumannOverlap_Plex(DM dm, IS *ovl, Mat *J, PetscErrorCode (**setup)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void **setup_ctx)
227*907a3e9cSStefano Zampini {
228*907a3e9cSStefano Zampini   DM                     odm;
229*907a3e9cSStefano Zampini   Mat                    pJ;
230*907a3e9cSStefano Zampini   PetscSF                sf = NULL;
231*907a3e9cSStefano Zampini   PetscSection           sec, osec;
232*907a3e9cSStefano Zampini   ISLocalToGlobalMapping l2g;
233*907a3e9cSStefano Zampini   const PetscInt        *idxs;
234*907a3e9cSStefano Zampini   PetscInt               n, mh;
235*907a3e9cSStefano Zampini 
236*907a3e9cSStefano Zampini   PetscFunctionBegin;
237*907a3e9cSStefano Zampini   *setup     = NULL;
238*907a3e9cSStefano Zampini   *setup_ctx = NULL;
239*907a3e9cSStefano Zampini   *ovl       = NULL;
240*907a3e9cSStefano Zampini   *J         = NULL;
241*907a3e9cSStefano Zampini 
242*907a3e9cSStefano Zampini   /* Overlapped mesh
243*907a3e9cSStefano Zampini      overlap is a little more generous, since it is not computed starting from the owned (Dirichlet) points, but from the locally owned cells */
244*907a3e9cSStefano Zampini   PetscCall(DMPlexDistributeOverlap(dm, 1, &sf, &odm));
245*907a3e9cSStefano Zampini   if (!odm) {
246*907a3e9cSStefano Zampini     PetscCall(PetscSFDestroy(&sf));
247*907a3e9cSStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
248*907a3e9cSStefano Zampini   }
249*907a3e9cSStefano Zampini 
250*907a3e9cSStefano Zampini   /* share discretization */
251*907a3e9cSStefano Zampini   PetscCall(DMGetLocalSection(dm, &sec));
252*907a3e9cSStefano Zampini   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)sec), &osec));
253*907a3e9cSStefano Zampini   PetscCall(PetscSFDistributeSection(sf, sec, NULL, osec));
254*907a3e9cSStefano Zampini   /* what to do here? using both is fine? */
255*907a3e9cSStefano Zampini   PetscCall(DMSetLocalSection(odm, osec));
256*907a3e9cSStefano Zampini   PetscCall(DMCopyDisc(dm, odm));
257*907a3e9cSStefano Zampini   PetscCall(DMPlexGetMaxProjectionHeight(dm, &mh));
258*907a3e9cSStefano Zampini   PetscCall(DMPlexSetMaxProjectionHeight(odm, mh));
259*907a3e9cSStefano Zampini   PetscCall(PetscSectionDestroy(&osec));
260*907a3e9cSStefano Zampini 
261*907a3e9cSStefano Zampini   /* material parameters */
262*907a3e9cSStefano Zampini   {
263*907a3e9cSStefano Zampini     Vec A;
264*907a3e9cSStefano Zampini 
265*907a3e9cSStefano Zampini     PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &A));
266*907a3e9cSStefano Zampini     if (A) {
267*907a3e9cSStefano Zampini       DM  dmAux, ocdm, odmAux;
268*907a3e9cSStefano Zampini       Vec oA;
269*907a3e9cSStefano Zampini 
270*907a3e9cSStefano Zampini       PetscCall(VecGetDM(A, &dmAux));
271*907a3e9cSStefano Zampini       PetscCall(DMClone(odm, &odmAux));
272*907a3e9cSStefano Zampini       PetscCall(DMGetCoordinateDM(odm, &ocdm));
273*907a3e9cSStefano Zampini       PetscCall(DMSetCoordinateDM(odmAux, ocdm));
274*907a3e9cSStefano Zampini       PetscCall(DMCopyDisc(dmAux, odmAux));
275*907a3e9cSStefano Zampini 
276*907a3e9cSStefano Zampini       PetscCall(DMGetLocalSection(dmAux, &sec));
277*907a3e9cSStefano Zampini       PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)sec), &osec));
278*907a3e9cSStefano Zampini       PetscCall(VecCreate(PetscObjectComm((PetscObject)A), &oA));
279*907a3e9cSStefano Zampini       PetscCall(VecSetDM(oA, odmAux));
280*907a3e9cSStefano Zampini       /* TODO: what if these values changes? */
281*907a3e9cSStefano Zampini       PetscCall(DMPlexDistributeField(dmAux, sf, sec, A, osec, oA));
282*907a3e9cSStefano Zampini       PetscCall(DMSetLocalSection(odmAux, osec));
283*907a3e9cSStefano Zampini       PetscCall(PetscSectionDestroy(&osec));
284*907a3e9cSStefano Zampini       PetscCall(DMSetAuxiliaryVec(odm, NULL, 0, 0, oA));
285*907a3e9cSStefano Zampini       PetscCall(VecDestroy(&oA));
286*907a3e9cSStefano Zampini       PetscCall(DMDestroy(&odmAux));
287*907a3e9cSStefano Zampini     }
288*907a3e9cSStefano Zampini   }
289*907a3e9cSStefano Zampini   PetscCall(PetscSFDestroy(&sf));
290*907a3e9cSStefano Zampini 
291*907a3e9cSStefano Zampini   PetscCall(DMViewFromOptions(dm, NULL, "-dm_plex_view_neumann_original"));
292*907a3e9cSStefano Zampini   PetscCall(PetscObjectSetName((PetscObject)odm, "OVL"));
293*907a3e9cSStefano Zampini   PetscCall(DMViewFromOptions(odm, NULL, "-dm_plex_view_neumann_overlap"));
294*907a3e9cSStefano Zampini 
295*907a3e9cSStefano Zampini   /* MATIS for the overlap region
296*907a3e9cSStefano Zampini      the HPDDM interface wants local matrices,
297*907a3e9cSStefano Zampini      we attach the global MATIS to the overlap IS,
298*907a3e9cSStefano Zampini      since we need it to do assembly */
299*907a3e9cSStefano Zampini   PetscCall(DMSetMatType(odm, MATIS));
300*907a3e9cSStefano Zampini   PetscCall(DMCreateMatrix(odm, &pJ));
301*907a3e9cSStefano Zampini   PetscCall(MatISGetLocalMat(pJ, J));
302*907a3e9cSStefano Zampini   PetscCall(PetscObjectReference((PetscObject)*J));
303*907a3e9cSStefano Zampini 
304*907a3e9cSStefano Zampini   /* overlap IS */
305*907a3e9cSStefano Zampini   PetscCall(MatISGetLocalToGlobalMapping(pJ, &l2g, NULL));
306*907a3e9cSStefano Zampini   PetscCall(ISLocalToGlobalMappingGetSize(l2g, &n));
307*907a3e9cSStefano Zampini   PetscCall(ISLocalToGlobalMappingGetIndices(l2g, &idxs));
308*907a3e9cSStefano Zampini   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)odm), n, idxs, PETSC_COPY_VALUES, ovl));
309*907a3e9cSStefano Zampini   PetscCall(ISLocalToGlobalMappingRestoreIndices(l2g, &idxs));
310*907a3e9cSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)*ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject)pJ));
311*907a3e9cSStefano Zampini   PetscCall(DMDestroy(&odm));
312*907a3e9cSStefano Zampini   PetscCall(MatDestroy(&pJ));
313*907a3e9cSStefano Zampini 
314*907a3e9cSStefano Zampini   /* special purpose setup function (composed in DMPlexSetSNESLocalFEM) */
315*907a3e9cSStefano Zampini   PetscCall(PetscObjectQueryFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", setup));
316*907a3e9cSStefano Zampini   if (*setup) PetscCall(PetscObjectCompose((PetscObject)*ovl, "_DM_Original_HPDDM", (PetscObject)dm));
317*907a3e9cSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
318*907a3e9cSStefano Zampini }
319