xref: /petsc/src/dm/impls/plex/plexpreallocate.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>
30c312b8eSJed Brown #include <petscsf.h>
4a9fb6443SMatthew G. Knepley #include <petscds.h>
5cb1e1211SMatthew G Knepley 
676185916SToby Isaac /* get adjacencies due to point-to-point constraints that can't be found with DMPlexGetAdjacency() */
79371c9d4SSatish Balay static PetscErrorCode DMPlexComputeAnchorAdjacencies(DM dm, PetscBool useCone, PetscBool useClosure, PetscSection *anchorSectionAdj, PetscInt *anchorAdj[]) {
876185916SToby Isaac   PetscInt     pStart, pEnd;
9e101f074SMatthew G. Knepley   PetscSection section, sectionGlobal, adjSec, aSec;
1076185916SToby Isaac   IS           aIS;
1176185916SToby Isaac 
1276185916SToby Isaac   PetscFunctionBegin;
139566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
149566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
159566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)section), &adjSec));
169566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
179566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(adjSec, pStart, pEnd));
1876185916SToby Isaac 
199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS));
2076185916SToby Isaac   if (aSec) {
2176185916SToby Isaac     const PetscInt *anchors;
2276185916SToby Isaac     PetscInt        p, q, a, aSize, *offsets, aStart, aEnd, *inverse, iSize, *adj, adjSize;
2376185916SToby Isaac     PetscInt       *tmpAdjP = NULL, *tmpAdjQ = NULL;
2476185916SToby Isaac     PetscSection    inverseSec;
2576185916SToby Isaac 
2676185916SToby Isaac     /* invert the constraint-to-anchor map */
279566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)aSec), &inverseSec));
289566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(inverseSec, pStart, pEnd));
299566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(aIS, &aSize));
309566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(aIS, &anchors));
3176185916SToby Isaac 
3276185916SToby Isaac     for (p = 0; p < aSize; p++) {
3376185916SToby Isaac       PetscInt a = anchors[p];
3476185916SToby Isaac 
359566063dSJacob Faibussowitsch       PetscCall(PetscSectionAddDof(inverseSec, a, 1));
3676185916SToby Isaac     }
379566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(inverseSec));
389566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(inverseSec, &iSize));
399566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(iSize, &inverse));
409566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(pEnd - pStart, &offsets));
419566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
4276185916SToby Isaac     for (p = aStart; p < aEnd; p++) {
4376185916SToby Isaac       PetscInt dof, off;
4476185916SToby Isaac 
459566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(aSec, p, &dof));
469566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(aSec, p, &off));
4776185916SToby Isaac 
4876185916SToby Isaac       for (q = 0; q < dof; q++) {
4976185916SToby Isaac         PetscInt iOff;
5076185916SToby Isaac 
5176185916SToby Isaac         a = anchors[off + q];
529566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(inverseSec, a, &iOff));
5376185916SToby Isaac         inverse[iOff + offsets[a - pStart]++] = p;
5476185916SToby Isaac       }
5576185916SToby Isaac     }
569566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(aIS, &anchors));
579566063dSJacob Faibussowitsch     PetscCall(PetscFree(offsets));
5876185916SToby Isaac 
5976185916SToby Isaac     /* construct anchorAdj and adjSec
6076185916SToby Isaac      *
6176185916SToby Isaac      * loop over anchors:
6276185916SToby Isaac      *   construct anchor adjacency
6376185916SToby Isaac      *   loop over constrained:
6476185916SToby Isaac      *     construct constrained adjacency
6576185916SToby Isaac      *     if not in anchor adjacency, add to dofs
6676185916SToby Isaac      * setup adjSec, allocate anchorAdj
6776185916SToby Isaac      * loop over anchors:
6876185916SToby Isaac      *   construct anchor adjacency
6976185916SToby Isaac      *   loop over constrained:
7076185916SToby Isaac      *     construct constrained adjacency
7176185916SToby Isaac      *     if not in anchor adjacency
7276185916SToby Isaac      *       if not already in list, put in list
7376185916SToby Isaac      *   sort, unique, reduce dof count
7476185916SToby Isaac      * optional: compactify
7576185916SToby Isaac      */
7676185916SToby Isaac     for (p = pStart; p < pEnd; p++) {
7776185916SToby Isaac       PetscInt iDof, iOff, i, r, s, numAdjP = PETSC_DETERMINE;
7876185916SToby Isaac 
799566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(inverseSec, p, &iDof));
8076185916SToby Isaac       if (!iDof) continue;
819566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(inverseSec, p, &iOff));
829566063dSJacob Faibussowitsch       PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, PETSC_TRUE, &numAdjP, &tmpAdjP));
8376185916SToby Isaac       for (i = 0; i < iDof; i++) {
8476185916SToby Isaac         PetscInt iNew = 0, qAdj, qAdjDof, qAdjCDof, numAdjQ = PETSC_DETERMINE;
8576185916SToby Isaac 
8676185916SToby Isaac         q = inverse[iOff + i];
879566063dSJacob Faibussowitsch         PetscCall(DMPlexGetAdjacency_Internal(dm, q, useCone, useClosure, PETSC_TRUE, &numAdjQ, &tmpAdjQ));
8876185916SToby Isaac         for (r = 0; r < numAdjQ; r++) {
8976185916SToby Isaac           qAdj = tmpAdjQ[r];
9076185916SToby Isaac           if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
9176185916SToby Isaac           for (s = 0; s < numAdjP; s++) {
9276185916SToby Isaac             if (qAdj == tmpAdjP[s]) break;
9376185916SToby Isaac           }
9476185916SToby Isaac           if (s < numAdjP) continue;
959566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(section, qAdj, &qAdjDof));
969566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetConstraintDof(section, qAdj, &qAdjCDof));
9776185916SToby Isaac           iNew += qAdjDof - qAdjCDof;
9876185916SToby Isaac         }
999566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(adjSec, p, iNew));
10076185916SToby Isaac       }
10176185916SToby Isaac     }
10276185916SToby Isaac 
1039566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(adjSec));
1049566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(adjSec, &adjSize));
1059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(adjSize, &adj));
10676185916SToby Isaac 
10776185916SToby Isaac     for (p = pStart; p < pEnd; p++) {
10876185916SToby Isaac       PetscInt iDof, iOff, i, r, s, aOff, aOffOrig, aDof, numAdjP = PETSC_DETERMINE;
10976185916SToby Isaac 
1109566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(inverseSec, p, &iDof));
11176185916SToby Isaac       if (!iDof) continue;
1129566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(inverseSec, p, &iOff));
1139566063dSJacob Faibussowitsch       PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, PETSC_TRUE, &numAdjP, &tmpAdjP));
1149566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(adjSec, p, &aDof));
1159566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(adjSec, p, &aOff));
11676185916SToby Isaac       aOffOrig = aOff;
11776185916SToby Isaac       for (i = 0; i < iDof; i++) {
11876185916SToby Isaac         PetscInt qAdj, qAdjDof, qAdjCDof, qAdjOff, nd, numAdjQ = PETSC_DETERMINE;
11976185916SToby Isaac 
12076185916SToby Isaac         q = inverse[iOff + i];
1219566063dSJacob Faibussowitsch         PetscCall(DMPlexGetAdjacency_Internal(dm, q, useCone, useClosure, PETSC_TRUE, &numAdjQ, &tmpAdjQ));
12276185916SToby Isaac         for (r = 0; r < numAdjQ; r++) {
12376185916SToby Isaac           qAdj = tmpAdjQ[r];
12476185916SToby Isaac           if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
12576185916SToby Isaac           for (s = 0; s < numAdjP; s++) {
12676185916SToby Isaac             if (qAdj == tmpAdjP[s]) break;
12776185916SToby Isaac           }
12876185916SToby Isaac           if (s < numAdjP) continue;
1299566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(section, qAdj, &qAdjDof));
1309566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetConstraintDof(section, qAdj, &qAdjCDof));
1319566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetOffset(sectionGlobal, qAdj, &qAdjOff));
1329371c9d4SSatish Balay           for (nd = 0; nd < qAdjDof - qAdjCDof; ++nd) { adj[aOff++] = (qAdjOff < 0 ? -(qAdjOff + 1) : qAdjOff) + nd; }
13376185916SToby Isaac         }
13476185916SToby Isaac       }
1359566063dSJacob Faibussowitsch       PetscCall(PetscSortRemoveDupsInt(&aDof, &adj[aOffOrig]));
1369566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(adjSec, p, aDof));
13776185916SToby Isaac     }
13876185916SToby Isaac     *anchorAdj = adj;
13976185916SToby Isaac 
14076185916SToby Isaac     /* clean up */
1419566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&inverseSec));
1429566063dSJacob Faibussowitsch     PetscCall(PetscFree(inverse));
1439566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmpAdjP));
1449566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmpAdjQ));
1459371c9d4SSatish Balay   } else {
14676185916SToby Isaac     *anchorAdj = NULL;
1479566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(adjSec));
14876185916SToby Isaac   }
14976185916SToby Isaac   *anchorSectionAdj = adjSec;
15076185916SToby Isaac   PetscFunctionReturn(0);
15176185916SToby Isaac }
15276185916SToby Isaac 
1539371c9d4SSatish Balay static PetscErrorCode DMPlexCreateAdjacencySection_Static(DM dm, PetscInt bs, PetscSF sfDof, PetscBool useCone, PetscBool useClosure, PetscBool useAnchors, PetscSection *sA, PetscInt **colIdx) {
154cb1e1211SMatthew G Knepley   MPI_Comm           comm;
155a9fb6443SMatthew G. Knepley   PetscMPIInt        size;
156a9fb6443SMatthew G. Knepley   PetscBool          doCommLocal, doComm, debug = PETSC_FALSE;
157a9fb6443SMatthew G. Knepley   PetscSF            sf, sfAdj;
158e101f074SMatthew G. Knepley   PetscSection       section, sectionGlobal, leafSectionAdj, rootSectionAdj, sectionAdj, anchorSectionAdj;
159a9fb6443SMatthew G. Knepley   PetscInt           nroots, nleaves, l, p, r;
160cb1e1211SMatthew G Knepley   const PetscInt    *leaves;
161cb1e1211SMatthew G Knepley   const PetscSFNode *remotes;
162cb1e1211SMatthew G Knepley   PetscInt           dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols;
1638821704fSMatthew G. Knepley   PetscInt          *tmpAdj = NULL, *adj, *rootAdj, *anchorAdj = NULL, *cols, *remoteOffsets;
16470034214SMatthew G. Knepley   PetscInt           adjSize;
165cb1e1211SMatthew G Knepley 
166cb1e1211SMatthew G Knepley   PetscFunctionBegin;
1679566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_view_preallocation", &debug, NULL));
1699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1709566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1719566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
1729566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
1739566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
1749566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
17547a6104aSMatthew G. Knepley   doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE;
1761c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm));
177cb1e1211SMatthew G Knepley   /* Create section for dof adjacency (dof ==> # adj dof) */
1789566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1799566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &numDof));
1809566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &leafSectionAdj));
1819566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(leafSectionAdj, 0, numDof));
1829566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &rootSectionAdj));
1839566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(rootSectionAdj, 0, numDof));
184cb1e1211SMatthew G Knepley   /*   Fill in the ghost dofs on the interface */
1859566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes));
186cb1e1211SMatthew G Knepley   /*
187964bf7afSMatthew G. Knepley    section        - maps points to (# dofs, local dofs)
188964bf7afSMatthew G. Knepley    sectionGlobal  - maps points to (# dofs, global dofs)
189964bf7afSMatthew G. Knepley    leafSectionAdj - maps unowned local dofs to # adj dofs
190964bf7afSMatthew G. Knepley    rootSectionAdj - maps   owned local dofs to # adj dofs
191964bf7afSMatthew G. Knepley    adj            - adj global dofs indexed by leafSectionAdj
192964bf7afSMatthew G. Knepley    rootAdj        - adj global dofs indexed by rootSectionAdj
193964bf7afSMatthew G. Knepley    sf    - describes shared points across procs
194964bf7afSMatthew G. Knepley    sfDof - describes shared dofs across procs
195964bf7afSMatthew G. Knepley    sfAdj - describes shared adjacent dofs across procs
196cb1e1211SMatthew G Knepley    ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point.
19776185916SToby Isaac   (0). If there are point-to-point constraints, add the adjacencies of constrained points to anchors in anchorAdj
19876185916SToby Isaac        (This is done in DMPlexComputeAnchorAdjacencies())
199cb1e1211SMatthew G Knepley     1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj
200cb1e1211SMatthew G Knepley        Reduce those counts to rootSectionAdj (now redundantly counting some interface points)
201cb1e1211SMatthew G Knepley     2. Visit owned points on interface, count adjacencies placing in rootSectionAdj
202cb1e1211SMatthew G Knepley        Create sfAdj connecting rootSectionAdj and leafSectionAdj
203cb1e1211SMatthew G Knepley     3. Visit unowned points on interface, write adjacencies to adj
204cb1e1211SMatthew G Knepley        Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies)
205cb1e1211SMatthew G Knepley     4. Visit owned points on interface, write adjacencies to rootAdj
206cb1e1211SMatthew G Knepley        Remove redundancy in rootAdj
207cb1e1211SMatthew G Knepley    ** The last two traversals use transitive closure
208cb1e1211SMatthew G Knepley     5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj)
209cb1e1211SMatthew G Knepley        Allocate memory addressed by sectionAdj (cols)
210cb1e1211SMatthew G Knepley     6. Visit all owned points in the subdomain, insert dof adjacencies into cols
211cb1e1211SMatthew G Knepley    ** Knowing all the column adjacencies, check ownership and sum into dnz and onz
212cb1e1211SMatthew G Knepley   */
2139566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeAnchorAdjacencies(dm, useCone, useClosure, &anchorSectionAdj, &anchorAdj));
214cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
21576185916SToby Isaac     PetscInt dof, off, d, q, anDof;
21670034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
217cb1e1211SMatthew G Knepley 
218fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
2199566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
2209566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
2219566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
222cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
223cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
224cb1e1211SMatthew G Knepley       PetscInt       ndof, ncdof;
225cb1e1211SMatthew G Knepley 
226cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
2279566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
2289566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
229*48a46eb9SPierre Jolivet       for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(leafSectionAdj, d, ndof - ncdof));
230cb1e1211SMatthew G Knepley     }
2319566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
23276185916SToby Isaac     if (anDof) {
233*48a46eb9SPierre Jolivet       for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(leafSectionAdj, d, anDof));
23476185916SToby Isaac     }
235cb1e1211SMatthew G Knepley   }
2369566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSectionAdj));
237cb1e1211SMatthew G Knepley   if (debug) {
2389566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n"));
2399566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(leafSectionAdj, NULL));
240cb1e1211SMatthew G Knepley   }
241cb1e1211SMatthew G Knepley   /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */
24247a6104aSMatthew G. Knepley   if (doComm) {
2439566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
2449566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
24569c11d05SVaclav Hapla     PetscCall(PetscSectionInvalidateMaxDof_Internal(rootSectionAdj));
246cb1e1211SMatthew G Knepley   }
247cb1e1211SMatthew G Knepley   if (debug) {
2489566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n"));
2499566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
250cb1e1211SMatthew G Knepley   }
251cb1e1211SMatthew G Knepley   /* Add in local adjacency sizes for owned dofs on interface (roots) */
252cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
25376185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof;
254cb1e1211SMatthew G Knepley 
2559566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
2569566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
257cb1e1211SMatthew G Knepley     if (!dof) continue;
2589566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
259cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
2609566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
261cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
262cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
263cb1e1211SMatthew G Knepley       PetscInt       ndof, ncdof;
264cb1e1211SMatthew G Knepley 
265cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
2669566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
2679566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
268*48a46eb9SPierre Jolivet       for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(rootSectionAdj, d, ndof - ncdof));
269cb1e1211SMatthew G Knepley     }
2709566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
27176185916SToby Isaac     if (anDof) {
272*48a46eb9SPierre Jolivet       for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(rootSectionAdj, d, anDof));
27376185916SToby Isaac     }
274cb1e1211SMatthew G Knepley   }
2759566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(rootSectionAdj));
276cb1e1211SMatthew G Knepley   if (debug) {
2779566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n"));
2789566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
279cb1e1211SMatthew G Knepley   }
280cb1e1211SMatthew G Knepley   /* Create adj SF based on dof SF */
2819566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets));
2829566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj));
2839566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
2846cded2eaSMatthew G. Knepley   if (debug && size > 1) {
2859566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency SF for Preallocation:\n"));
2869566063dSJacob Faibussowitsch     PetscCall(PetscSFView(sfAdj, NULL));
287cb1e1211SMatthew G Knepley   }
288cb1e1211SMatthew G Knepley   /* Create leaf adjacency */
2899566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSectionAdj));
2909566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(leafSectionAdj, &adjSize));
2919566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(adjSize, &adj));
292cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
29376185916SToby Isaac     PetscInt dof, off, d, q, anDof, anOff;
29470034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
295cb1e1211SMatthew G Knepley 
296fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
2979566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
2989566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
2999566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
3009566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
3019566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
302cb1e1211SMatthew G Knepley     for (d = off; d < off + dof; ++d) {
303cb1e1211SMatthew G Knepley       PetscInt aoff, i = 0;
304cb1e1211SMatthew G Knepley 
3059566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(leafSectionAdj, d, &aoff));
306cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
307cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
308cb1e1211SMatthew G Knepley         PetscInt       ndof, ncdof, ngoff, nd;
309cb1e1211SMatthew G Knepley 
310cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
3119566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
3129566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
3139566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
314cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof - ncdof; ++nd) {
315cb1e1211SMatthew G Knepley           adj[aoff + i] = (ngoff < 0 ? -(ngoff + 1) : ngoff) + nd;
316cb1e1211SMatthew G Knepley           ++i;
317cb1e1211SMatthew G Knepley         }
318cb1e1211SMatthew G Knepley       }
31976185916SToby Isaac       for (q = 0; q < anDof; q++) {
32076185916SToby Isaac         adj[aoff + i] = anchorAdj[anOff + q];
32176185916SToby Isaac         ++i;
32276185916SToby Isaac       }
323cb1e1211SMatthew G Knepley     }
324cb1e1211SMatthew G Knepley   }
325cb1e1211SMatthew G Knepley   /* Debugging */
326cb1e1211SMatthew G Knepley   if (debug) {
327cb1e1211SMatthew G Knepley     IS tmp;
3289566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Leaf adjacency indices\n"));
3299566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp));
3309566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
3319566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
332cb1e1211SMatthew G Knepley   }
333543482b8SMatthew G. Knepley   /* Gather adjacent indices to root */
3349566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(rootSectionAdj, &adjSize));
3359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(adjSize, &rootAdj));
336cb1e1211SMatthew G Knepley   for (r = 0; r < adjSize; ++r) rootAdj[r] = -1;
33747a6104aSMatthew G. Knepley   if (doComm) {
338543482b8SMatthew G. Knepley     const PetscInt *indegree;
339543482b8SMatthew G. Knepley     PetscInt       *remoteadj, radjsize = 0;
340543482b8SMatthew G. Knepley 
3419566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sfAdj, &indegree));
3429566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sfAdj, &indegree));
343543482b8SMatthew G. Knepley     for (p = 0; p < adjSize; ++p) radjsize += indegree[p];
3449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(radjsize, &remoteadj));
3459566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj));
3469566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherEnd(sfAdj, MPIU_INT, adj, remoteadj));
3476dba2905SMatthew G. Knepley     for (p = 0, l = 0, r = 0; p < adjSize; ++p, l = PetscMax(p, l + indegree[p - 1])) {
348543482b8SMatthew G. Knepley       PetscInt s;
3496dba2905SMatthew G. Knepley       for (s = 0; s < indegree[p]; ++s, ++r) rootAdj[l + s] = remoteadj[r];
350543482b8SMatthew G. Knepley     }
35163a3b9bcSJacob Faibussowitsch     PetscCheck(r == radjsize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %" PetscInt_FMT " != %" PetscInt_FMT, r, radjsize);
35263a3b9bcSJacob Faibussowitsch     PetscCheck(l == adjSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %" PetscInt_FMT " != %" PetscInt_FMT, l, adjSize);
3539566063dSJacob Faibussowitsch     PetscCall(PetscFree(remoteadj));
354cb1e1211SMatthew G Knepley   }
3559566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfAdj));
3569566063dSJacob Faibussowitsch   PetscCall(PetscFree(adj));
357cb1e1211SMatthew G Knepley   /* Debugging */
358cb1e1211SMatthew G Knepley   if (debug) {
359cb1e1211SMatthew G Knepley     IS tmp;
3609566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices after gather\n"));
3619566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
3629566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
3639566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
364cb1e1211SMatthew G Knepley   }
365cb1e1211SMatthew G Knepley   /* Add in local adjacency indices for owned dofs on interface (roots) */
366cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
36776185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff;
368cb1e1211SMatthew G Knepley 
3699566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
3709566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
371cb1e1211SMatthew G Knepley     if (!dof) continue;
3729566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
373cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
3749566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
3759566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
3769566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
377cb1e1211SMatthew G Knepley     for (d = off; d < off + dof; ++d) {
378cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i;
379cb1e1211SMatthew G Knepley 
3809566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof));
3819566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
382cb1e1211SMatthew G Knepley       i = adof - 1;
38376185916SToby Isaac       for (q = 0; q < anDof; q++) {
38476185916SToby Isaac         rootAdj[aoff + i] = anchorAdj[anOff + q];
38576185916SToby Isaac         --i;
38676185916SToby Isaac       }
387cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
388cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
389cb1e1211SMatthew G Knepley         PetscInt       ndof, ncdof, ngoff, nd;
390cb1e1211SMatthew G Knepley 
391cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
3929566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
3939566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
3949566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
395cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof - ncdof; ++nd) {
396cb1e1211SMatthew G Knepley           rootAdj[aoff + i] = ngoff < 0 ? -(ngoff + 1) + nd : ngoff + nd;
397cb1e1211SMatthew G Knepley           --i;
398cb1e1211SMatthew G Knepley         }
399cb1e1211SMatthew G Knepley       }
400cb1e1211SMatthew G Knepley     }
401cb1e1211SMatthew G Knepley   }
402cb1e1211SMatthew G Knepley   /* Debugging */
403cb1e1211SMatthew G Knepley   if (debug) {
404cb1e1211SMatthew G Knepley     IS tmp;
4059566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices\n"));
4069566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
4079566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
4089566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
409cb1e1211SMatthew G Knepley   }
410cb1e1211SMatthew G Knepley   /* Compress indices */
4119566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(rootSectionAdj));
412cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
413cb1e1211SMatthew G Knepley     PetscInt dof, cdof, off, d;
414cb1e1211SMatthew G Knepley     PetscInt adof, aoff;
415cb1e1211SMatthew G Knepley 
4169566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4179566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
4189566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
419cb1e1211SMatthew G Knepley     if (!dof) continue;
4209566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
421cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
422cb1e1211SMatthew G Knepley     for (d = off; d < off + dof - cdof; ++d) {
4239566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof));
4249566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
4259566063dSJacob Faibussowitsch       PetscCall(PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]));
4269566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(rootSectionAdj, d, adof));
427cb1e1211SMatthew G Knepley     }
428cb1e1211SMatthew G Knepley   }
429cb1e1211SMatthew G Knepley   /* Debugging */
430cb1e1211SMatthew G Knepley   if (debug) {
431cb1e1211SMatthew G Knepley     IS tmp;
4329566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n"));
4339566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
4349566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices after compression\n"));
4359566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
4369566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
4379566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
438cb1e1211SMatthew G Knepley   }
439cb1e1211SMatthew G Knepley   /* Build adjacency section: Maps global indices to sets of adjacent global indices */
4409566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd));
4419566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &sectionAdj));
4429566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd));
443cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
44476185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof;
445cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
446cb1e1211SMatthew G Knepley 
4479566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4489566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
4499566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
4509566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
451cb1e1211SMatthew G Knepley     for (d = 0; d < dof - cdof; ++d) {
452cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
453cb1e1211SMatthew G Knepley 
4549566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSectionAdj, off + d, &ldof));
4559566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, off + d, &rdof));
456cb1e1211SMatthew G Knepley       if (ldof > 0) {
457cb1e1211SMatthew G Knepley         /* We do not own this point */
458cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
4599566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetDof(sectionAdj, goff + d, rdof));
460cb1e1211SMatthew G Knepley       } else {
461cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
462cb1e1211SMatthew G Knepley       }
463cb1e1211SMatthew G Knepley     }
464cb1e1211SMatthew G Knepley     if (found) continue;
4659566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4669566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
4679566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
468cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
469cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
470cb1e1211SMatthew G Knepley       PetscInt       ndof, ncdof, noff;
471cb1e1211SMatthew G Knepley 
472cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
4739566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
4749566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
4759566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, padj, &noff));
476*48a46eb9SPierre Jolivet       for (d = goff; d < goff + dof - cdof; ++d) PetscCall(PetscSectionAddDof(sectionAdj, d, ndof - ncdof));
477cb1e1211SMatthew G Knepley     }
4789566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
47976185916SToby Isaac     if (anDof) {
480*48a46eb9SPierre Jolivet       for (d = goff; d < goff + dof - cdof; ++d) PetscCall(PetscSectionAddDof(sectionAdj, d, anDof));
48176185916SToby Isaac     }
482cb1e1211SMatthew G Knepley   }
4839566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(sectionAdj));
484cb1e1211SMatthew G Knepley   if (debug) {
4859566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation:\n"));
4869566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(sectionAdj, NULL));
487cb1e1211SMatthew G Knepley   }
488cb1e1211SMatthew G Knepley   /* Get adjacent indices */
4899566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(sectionAdj, &numCols));
4909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numCols, &cols));
491cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
49276185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff;
493cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
494cb1e1211SMatthew G Knepley 
4959566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4969566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
4979566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
4989566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
499cb1e1211SMatthew G Knepley     for (d = 0; d < dof - cdof; ++d) {
500cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
501cb1e1211SMatthew G Knepley 
5029566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSectionAdj, off + d, &ldof));
5039566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, off + d, &rdof));
504cb1e1211SMatthew G Knepley       if (ldof > 0) {
505cb1e1211SMatthew G Knepley         /* We do not own this point */
506cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
507cb1e1211SMatthew G Knepley         PetscInt aoff, roff;
508cb1e1211SMatthew G Knepley 
5099566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, goff + d, &aoff));
5109566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(rootSectionAdj, off + d, &roff));
5119566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&cols[aoff], &rootAdj[roff], rdof));
512cb1e1211SMatthew G Knepley       } else {
513cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
514cb1e1211SMatthew G Knepley       }
515cb1e1211SMatthew G Knepley     }
516cb1e1211SMatthew G Knepley     if (found) continue;
5179566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
5189566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
5199566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
520cb1e1211SMatthew G Knepley     for (d = goff; d < goff + dof - cdof; ++d) {
521cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i = 0;
522cb1e1211SMatthew G Knepley 
5239566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, d, &adof));
5249566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(sectionAdj, d, &aoff));
525cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
526cb1e1211SMatthew G Knepley         const PetscInt  padj = tmpAdj[q];
527cb1e1211SMatthew G Knepley         PetscInt        ndof, ncdof, ngoff, nd;
528cb1e1211SMatthew G Knepley         const PetscInt *ncind;
529cb1e1211SMatthew G Knepley 
530cb1e1211SMatthew G Knepley         /* Adjacent points may not be in the section chart */
531cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
5329566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
5339566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
5349566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintIndices(section, padj, &ncind));
5359566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
5369371c9d4SSatish Balay         for (nd = 0; nd < ndof - ncdof; ++nd, ++i) { cols[aoff + i] = ngoff < 0 ? -(ngoff + 1) + nd : ngoff + nd; }
537cb1e1211SMatthew G Knepley       }
5389371c9d4SSatish Balay       for (q = 0; q < anDof; q++, i++) { cols[aoff + i] = anchorAdj[anOff + q]; }
53963a3b9bcSJacob Faibussowitsch       PetscCheck(i == adof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of entries %" PetscInt_FMT " != %" PetscInt_FMT " for dof %" PetscInt_FMT " (point %" PetscInt_FMT ")", i, adof, d, p);
540cb1e1211SMatthew G Knepley     }
541cb1e1211SMatthew G Knepley   }
5429566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&anchorSectionAdj));
5439566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&leafSectionAdj));
5449566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&rootSectionAdj));
5459566063dSJacob Faibussowitsch   PetscCall(PetscFree(anchorAdj));
5469566063dSJacob Faibussowitsch   PetscCall(PetscFree(rootAdj));
5479566063dSJacob Faibussowitsch   PetscCall(PetscFree(tmpAdj));
548cb1e1211SMatthew G Knepley   /* Debugging */
549cb1e1211SMatthew G Knepley   if (debug) {
550cb1e1211SMatthew G Knepley     IS tmp;
5519566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Column indices\n"));
5529566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp));
5539566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
5549566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
555cb1e1211SMatthew G Knepley   }
556a9fb6443SMatthew G. Knepley 
557a9fb6443SMatthew G. Knepley   *sA     = sectionAdj;
558a9fb6443SMatthew G. Knepley   *colIdx = cols;
559a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
560a9fb6443SMatthew G. Knepley }
561a9fb6443SMatthew G. Knepley 
5629371c9d4SSatish Balay static PetscErrorCode DMPlexUpdateAllocation_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[]) {
563e101f074SMatthew G. Knepley   PetscSection section;
564a9fb6443SMatthew G. Knepley   PetscInt     rStart, rEnd, r, pStart, pEnd, p;
565a9fb6443SMatthew G. Knepley 
566a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
567a9fb6443SMatthew G. Knepley   /* This loop needs to change to a loop over points, then field dofs, which means we need to look both sections */
5689566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
5691dca8a05SBarry Smith   PetscCheck((rStart % bs) == 0 && (rEnd % bs) == 0, PetscObjectComm((PetscObject)rLayout), PETSC_ERR_ARG_WRONG, "Invalid layout [%" PetscInt_FMT ", %" PetscInt_FMT ") for matrix, must be divisible by block size %" PetscInt_FMT, rStart, rEnd, bs);
570a9fb6443SMatthew G. Knepley   if (f >= 0 && bs == 1) {
5719566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
5729566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
573a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
574294b7009SMatthew G. Knepley       PetscInt rS, rE;
575a9fb6443SMatthew G. Knepley 
5769566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE));
577a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
578a9fb6443SMatthew G. Knepley         PetscInt numCols, cStart, c;
579a9fb6443SMatthew G. Knepley 
5809566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
5819566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
582a9fb6443SMatthew G. Knepley         for (c = cStart; c < cStart + numCols; ++c) {
583a9fb6443SMatthew G. Knepley           if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
584a9fb6443SMatthew G. Knepley             ++dnz[r - rStart];
585a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++dnzu[r - rStart];
586a9fb6443SMatthew G. Knepley           } else {
587a9fb6443SMatthew G. Knepley             ++onz[r - rStart];
588a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++onzu[r - rStart];
589a9fb6443SMatthew G. Knepley           }
590a9fb6443SMatthew G. Knepley         }
591a9fb6443SMatthew G. Knepley       }
592a9fb6443SMatthew G. Knepley     }
593a9fb6443SMatthew G. Knepley   } else {
594cb1e1211SMatthew G Knepley     /* Only loop over blocks of rows */
595cb1e1211SMatthew G Knepley     for (r = rStart / bs; r < rEnd / bs; ++r) {
596cb1e1211SMatthew G Knepley       const PetscInt row = r * bs;
597cb1e1211SMatthew G Knepley       PetscInt       numCols, cStart, c;
598cb1e1211SMatthew G Knepley 
5999566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, row, &numCols));
6009566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(sectionAdj, row, &cStart));
601cb1e1211SMatthew G Knepley       for (c = cStart; c < cStart + numCols; ++c) {
602e7bcfa36SMatthew G. Knepley         if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
603e7bcfa36SMatthew G. Knepley           ++dnz[r - rStart / bs];
604e7bcfa36SMatthew G. Knepley           if (cols[c] >= row) ++dnzu[r - rStart / bs];
605cb1e1211SMatthew G Knepley         } else {
606e7bcfa36SMatthew G. Knepley           ++onz[r - rStart / bs];
607e7bcfa36SMatthew G. Knepley           if (cols[c] >= row) ++onzu[r - rStart / bs];
608cb1e1211SMatthew G Knepley         }
609cb1e1211SMatthew G Knepley       }
610cb1e1211SMatthew G Knepley     }
611a9fb6443SMatthew G. Knepley     for (r = 0; r < (rEnd - rStart) / bs; ++r) {
612cb1e1211SMatthew G Knepley       dnz[r] /= bs;
613cb1e1211SMatthew G Knepley       onz[r] /= bs;
614cb1e1211SMatthew G Knepley       dnzu[r] /= bs;
615cb1e1211SMatthew G Knepley       onzu[r] /= bs;
616cb1e1211SMatthew G Knepley     }
617cb1e1211SMatthew G Knepley   }
618a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
619a9fb6443SMatthew G. Knepley }
620a9fb6443SMatthew G. Knepley 
6219371c9d4SSatish Balay static PetscErrorCode DMPlexFillMatrix_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], Mat A) {
622e101f074SMatthew G. Knepley   PetscSection section;
623a9fb6443SMatthew G. Knepley   PetscScalar *values;
624a9fb6443SMatthew G. Knepley   PetscInt     rStart, rEnd, r, pStart, pEnd, p, len, maxRowLen = 0;
625a9fb6443SMatthew G. Knepley 
626a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
6279566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
628a9fb6443SMatthew G. Knepley   for (r = rStart; r < rEnd; ++r) {
6299566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(sectionAdj, r, &len));
630a9fb6443SMatthew G. Knepley     maxRowLen = PetscMax(maxRowLen, len);
631a9fb6443SMatthew G. Knepley   }
6329566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(maxRowLen, &values));
633a9fb6443SMatthew G. Knepley   if (f >= 0 && bs == 1) {
6349566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
6359566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
636a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
637294b7009SMatthew G. Knepley       PetscInt rS, rE;
638a9fb6443SMatthew G. Knepley 
6399566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE));
640a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
6417e01ee02SMatthew G. Knepley         PetscInt numCols, cStart;
642a9fb6443SMatthew G. Knepley 
6439566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
6449566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
6459566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
646a9fb6443SMatthew G. Knepley       }
647a9fb6443SMatthew G. Knepley     }
648a9fb6443SMatthew G. Knepley   } else {
649a9fb6443SMatthew G. Knepley     for (r = rStart; r < rEnd; ++r) {
650a9fb6443SMatthew G. Knepley       PetscInt numCols, cStart;
651a9fb6443SMatthew G. Knepley 
6529566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
6539566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
6549566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
655a9fb6443SMatthew G. Knepley     }
656a9fb6443SMatthew G. Knepley   }
6579566063dSJacob Faibussowitsch   PetscCall(PetscFree(values));
658a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
659a9fb6443SMatthew G. Knepley }
660a9fb6443SMatthew G. Knepley 
66125e402d2SMatthew G. Knepley /*@C
66225e402d2SMatthew G. Knepley   DMPlexPreallocateOperator - Calculate the matrix nonzero pattern based upon the information in the DM,
66325e402d2SMatthew G. Knepley   the PetscDS it contains, and the default PetscSection.
66425e402d2SMatthew G. Knepley 
66525e402d2SMatthew G. Knepley   Collective
66625e402d2SMatthew G. Knepley 
6674165533cSJose E. Roman   Input Parameters:
66825e402d2SMatthew G. Knepley + dm   - The DMPlex
66925e402d2SMatthew G. Knepley . bs   - The matrix blocksize
67025e402d2SMatthew G. Knepley . dnz  - An array to hold the number of nonzeros in the diagonal block
67125e402d2SMatthew G. Knepley . onz  - An array to hold the number of nonzeros in the off-diagonal block
67225e402d2SMatthew G. Knepley . dnzu - An array to hold the number of nonzeros in the upper triangle of the diagonal block
67325e402d2SMatthew G. Knepley . onzu - An array to hold the number of nonzeros in the upper triangle of the off-diagonal block
674a2b725a8SWilliam Gropp - fillMatrix - If PETSC_TRUE, fill the matrix with zeros
67525e402d2SMatthew G. Knepley 
6764165533cSJose E. Roman   Output Parameter:
67725e402d2SMatthew G. Knepley . A - The preallocated matrix
67825e402d2SMatthew G. Knepley 
67925e402d2SMatthew G. Knepley   Level: advanced
68025e402d2SMatthew G. Knepley 
681db781477SPatrick Sanan .seealso: `DMCreateMatrix()`
68225e402d2SMatthew G. Knepley @*/
6839371c9d4SSatish Balay PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix) {
684a9fb6443SMatthew G. Knepley   MPI_Comm     comm;
685a9fb6443SMatthew G. Knepley   PetscDS      prob;
686a9fb6443SMatthew G. Knepley   MatType      mtype;
687a9fb6443SMatthew G. Knepley   PetscSF      sf, sfDof;
688e101f074SMatthew G. Knepley   PetscSection section;
689a9fb6443SMatthew G. Knepley   PetscInt    *remoteOffsets;
690a9fb6443SMatthew G. Knepley   PetscSection sectionAdj[4] = {NULL, NULL, NULL, NULL};
691a9fb6443SMatthew G. Knepley   PetscInt    *cols[4]       = {NULL, NULL, NULL, NULL};
692a9fb6443SMatthew G. Knepley   PetscBool    useCone, useClosure;
6937e01ee02SMatthew G. Knepley   PetscInt     Nf, f, idx, locRows;
694a9fb6443SMatthew G. Knepley   PetscLayout  rLayout;
695a9fb6443SMatthew G. Knepley   PetscBool    isSymBlock, isSymSeqBlock, isSymMPIBlock, debug = PETSC_FALSE;
6966cded2eaSMatthew G. Knepley   PetscMPIInt  size;
697a9fb6443SMatthew G. Knepley 
698a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
699a9fb6443SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
700064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(A, MAT_CLASSID, 7);
701dadcf809SJacob Faibussowitsch   if (dnz) PetscValidIntPointer(dnz, 3);
702dadcf809SJacob Faibussowitsch   if (onz) PetscValidIntPointer(onz, 4);
703dadcf809SJacob Faibussowitsch   if (dnzu) PetscValidIntPointer(dnzu, 5);
704dadcf809SJacob Faibussowitsch   if (onzu) PetscValidIntPointer(onzu, 6);
7059566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
7069566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
7079566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
7089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_view_preallocation", &debug, NULL));
7099566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
7109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
7119566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_Preallocate, dm, 0, 0, 0));
712a9fb6443SMatthew G. Knepley   /* Create dof SF based on point SF */
713a9fb6443SMatthew G. Knepley   if (debug) {
714e101f074SMatthew G. Knepley     PetscSection section, sectionGlobal;
715a9fb6443SMatthew G. Knepley     PetscSF      sf;
716a9fb6443SMatthew G. Knepley 
7179566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
7189566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
7199566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
7209566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Input Section for Preallocation:\n"));
7219566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(section, NULL));
7229566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Input Global Section for Preallocation:\n"));
7239566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(sectionGlobal, NULL));
7246cded2eaSMatthew G. Knepley     if (size > 1) {
7259566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Input SF for Preallocation:\n"));
7269566063dSJacob Faibussowitsch       PetscCall(PetscSFView(sf, NULL));
727a9fb6443SMatthew G. Knepley     }
7286cded2eaSMatthew G. Knepley   }
7299566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets));
7309566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof));
7319566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
7326cded2eaSMatthew G. Knepley   if (debug && size > 1) {
7339566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Dof SF for Preallocation:\n"));
7349566063dSJacob Faibussowitsch     PetscCall(PetscSFView(sfDof, NULL));
735a9fb6443SMatthew G. Knepley   }
736a9fb6443SMatthew G. Knepley   /* Create allocation vectors from adjacency graph */
7379566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &locRows, NULL));
7389566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
7399566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
7409566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
7419566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
742a9fb6443SMatthew G. Knepley   /* There are 4 types of adjacency */
7439566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
744a9fb6443SMatthew G. Knepley   if (Nf < 1 || bs > 1) {
7459566063dSJacob Faibussowitsch     PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
746a9fb6443SMatthew G. Knepley     idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7479566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]));
7489566063dSJacob Faibussowitsch     PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
749a9fb6443SMatthew G. Knepley   } else {
750a9fb6443SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
7519566063dSJacob Faibussowitsch       PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
752a9fb6443SMatthew G. Knepley       idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7539566063dSJacob Faibussowitsch       if (!sectionAdj[idx]) PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]));
7549566063dSJacob Faibussowitsch       PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
755a9fb6443SMatthew G. Knepley     }
756a9fb6443SMatthew G. Knepley   }
7579566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfDof));
758cb1e1211SMatthew G Knepley   /* Set matrix pattern */
7599566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu));
7609566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
76189545effSMatthew G. Knepley   /* Check for symmetric storage */
7629566063dSJacob Faibussowitsch   PetscCall(MatGetType(A, &mtype));
7639566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATSBAIJ, &isSymBlock));
7649566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock));
7659566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock));
7669566063dSJacob Faibussowitsch   if (isSymBlock || isSymSeqBlock || isSymMPIBlock) PetscCall(MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
767cb1e1211SMatthew G Knepley   /* Fill matrix with zeros */
768cb1e1211SMatthew G Knepley   if (fillMatrix) {
769a9fb6443SMatthew G. Knepley     if (Nf < 1 || bs > 1) {
7709566063dSJacob Faibussowitsch       PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
771a9fb6443SMatthew G. Knepley       idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7729566063dSJacob Faibussowitsch       PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], A));
773a9fb6443SMatthew G. Knepley     } else {
774a9fb6443SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
7759566063dSJacob Faibussowitsch         PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
776a9fb6443SMatthew G. Knepley         idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7779566063dSJacob Faibussowitsch         PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], A));
778cb1e1211SMatthew G Knepley       }
779cb1e1211SMatthew G Knepley     }
7809566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
7819566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
782cb1e1211SMatthew G Knepley   }
7839566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
7849371c9d4SSatish Balay   for (idx = 0; idx < 4; ++idx) {
7859371c9d4SSatish Balay     PetscCall(PetscSectionDestroy(&sectionAdj[idx]));
7869371c9d4SSatish Balay     PetscCall(PetscFree(cols[idx]));
7879371c9d4SSatish Balay   }
7889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_Preallocate, dm, 0, 0, 0));
789cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
790cb1e1211SMatthew G Knepley }
791cb1e1211SMatthew G Knepley 
792cb1e1211SMatthew G Knepley #if 0
793cb1e1211SMatthew G Knepley PetscErrorCode DMPlexPreallocateOperator_2(DM dm, PetscInt bs, PetscSection section, PetscSection sectionGlobal, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix)
794cb1e1211SMatthew G Knepley {
795cb1e1211SMatthew G Knepley   PetscInt       *tmpClosure,*tmpAdj,*visits;
796cb1e1211SMatthew G Knepley   PetscInt        c,cStart,cEnd,pStart,pEnd;
797cb1e1211SMatthew G Knepley 
798cb1e1211SMatthew G Knepley   PetscFunctionBegin;
7999566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
8009566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
8019566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
802cb1e1211SMatthew G Knepley 
803e7b80042SMatthew G. Knepley   maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1));
804cb1e1211SMatthew G Knepley 
8059566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
806cb1e1211SMatthew G Knepley   npoints = pEnd - pStart;
807cb1e1211SMatthew G Knepley 
8089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits));
8099566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(lvisits,pEnd-pStart));
8109566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(visits,pEnd-pStart));
8119566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
812cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
813cb1e1211SMatthew G Knepley     PetscInt *support = tmpClosure;
8149566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support));
815cb1e1211SMatthew G Knepley     for (p=0; p<supportSize; p++) lvisits[support[p]]++;
816cb1e1211SMatthew G Knepley   }
8179566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM));
8189566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd  (sf,MPIU_INT,lvisits,visits,MPI_SUM));
8199566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits,MPI_REPLACE));
8209566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd  (sf,MPIU_INT,visits,lvisits));
821cb1e1211SMatthew G Knepley 
8229566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks());
823cb1e1211SMatthew G Knepley 
8249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner));
825cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
8269566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(cellmat,maxClosureSize*maxClosureSize));
827cb1e1211SMatthew G Knepley     /*
828cb1e1211SMatthew G Knepley      Depth-first walk of transitive closure.
829cb1e1211SMatthew G Knepley      At each leaf frame f of transitive closure that we see, add 1/visits[f] to each pair (p,q) not marked as done in cellmat.
830cb1e1211SMatthew G Knepley      This contribution is added to dnz if owning ranks of p and q match, to onz otherwise.
831cb1e1211SMatthew G Knepley      */
832cb1e1211SMatthew G Knepley   }
833cb1e1211SMatthew G Knepley 
8349566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM));
8359566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd  (sf,MPIU_INT,lonz,onz,MPI_SUM));
836cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
837cb1e1211SMatthew G Knepley }
838cb1e1211SMatthew G Knepley #endif
839