xref: /petsc/src/dm/impls/plex/plexpreallocate.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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() */
7*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexComputeAnchorAdjacencies(DM dm, PetscBool useCone, PetscBool useClosure, PetscSection *anchorSectionAdj, PetscInt *anchorAdj[])
8*d71ae5a4SJacob Faibussowitsch {
976185916SToby Isaac   PetscInt     pStart, pEnd;
10e101f074SMatthew G. Knepley   PetscSection section, sectionGlobal, adjSec, aSec;
1176185916SToby Isaac   IS           aIS;
1276185916SToby Isaac 
1376185916SToby Isaac   PetscFunctionBegin;
149566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
159566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
169566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)section), &adjSec));
179566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
189566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(adjSec, pStart, pEnd));
1976185916SToby Isaac 
209566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS));
2176185916SToby Isaac   if (aSec) {
2276185916SToby Isaac     const PetscInt *anchors;
2376185916SToby Isaac     PetscInt        p, q, a, aSize, *offsets, aStart, aEnd, *inverse, iSize, *adj, adjSize;
2476185916SToby Isaac     PetscInt       *tmpAdjP = NULL, *tmpAdjQ = NULL;
2576185916SToby Isaac     PetscSection    inverseSec;
2676185916SToby Isaac 
2776185916SToby Isaac     /* invert the constraint-to-anchor map */
289566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)aSec), &inverseSec));
299566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(inverseSec, pStart, pEnd));
309566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(aIS, &aSize));
319566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(aIS, &anchors));
3276185916SToby Isaac 
3376185916SToby Isaac     for (p = 0; p < aSize; p++) {
3476185916SToby Isaac       PetscInt a = anchors[p];
3576185916SToby Isaac 
369566063dSJacob Faibussowitsch       PetscCall(PetscSectionAddDof(inverseSec, a, 1));
3776185916SToby Isaac     }
389566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(inverseSec));
399566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(inverseSec, &iSize));
409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(iSize, &inverse));
419566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(pEnd - pStart, &offsets));
429566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
4376185916SToby Isaac     for (p = aStart; p < aEnd; p++) {
4476185916SToby Isaac       PetscInt dof, off;
4576185916SToby Isaac 
469566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(aSec, p, &dof));
479566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(aSec, p, &off));
4876185916SToby Isaac 
4976185916SToby Isaac       for (q = 0; q < dof; q++) {
5076185916SToby Isaac         PetscInt iOff;
5176185916SToby Isaac 
5276185916SToby Isaac         a = anchors[off + q];
539566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(inverseSec, a, &iOff));
5476185916SToby Isaac         inverse[iOff + offsets[a - pStart]++] = p;
5576185916SToby Isaac       }
5676185916SToby Isaac     }
579566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(aIS, &anchors));
589566063dSJacob Faibussowitsch     PetscCall(PetscFree(offsets));
5976185916SToby Isaac 
6076185916SToby Isaac     /* construct anchorAdj and adjSec
6176185916SToby Isaac      *
6276185916SToby Isaac      * loop over anchors:
6376185916SToby Isaac      *   construct anchor adjacency
6476185916SToby Isaac      *   loop over constrained:
6576185916SToby Isaac      *     construct constrained adjacency
6676185916SToby Isaac      *     if not in anchor adjacency, add to dofs
6776185916SToby Isaac      * setup adjSec, allocate anchorAdj
6876185916SToby Isaac      * loop over anchors:
6976185916SToby Isaac      *   construct anchor adjacency
7076185916SToby Isaac      *   loop over constrained:
7176185916SToby Isaac      *     construct constrained adjacency
7276185916SToby Isaac      *     if not in anchor adjacency
7376185916SToby Isaac      *       if not already in list, put in list
7476185916SToby Isaac      *   sort, unique, reduce dof count
7576185916SToby Isaac      * optional: compactify
7676185916SToby Isaac      */
7776185916SToby Isaac     for (p = pStart; p < pEnd; p++) {
7876185916SToby Isaac       PetscInt iDof, iOff, i, r, s, numAdjP = PETSC_DETERMINE;
7976185916SToby Isaac 
809566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(inverseSec, p, &iDof));
8176185916SToby Isaac       if (!iDof) continue;
829566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(inverseSec, p, &iOff));
839566063dSJacob Faibussowitsch       PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, PETSC_TRUE, &numAdjP, &tmpAdjP));
8476185916SToby Isaac       for (i = 0; i < iDof; i++) {
8576185916SToby Isaac         PetscInt iNew = 0, qAdj, qAdjDof, qAdjCDof, numAdjQ = PETSC_DETERMINE;
8676185916SToby Isaac 
8776185916SToby Isaac         q = inverse[iOff + i];
889566063dSJacob Faibussowitsch         PetscCall(DMPlexGetAdjacency_Internal(dm, q, useCone, useClosure, PETSC_TRUE, &numAdjQ, &tmpAdjQ));
8976185916SToby Isaac         for (r = 0; r < numAdjQ; r++) {
9076185916SToby Isaac           qAdj = tmpAdjQ[r];
9176185916SToby Isaac           if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
9276185916SToby Isaac           for (s = 0; s < numAdjP; s++) {
9376185916SToby Isaac             if (qAdj == tmpAdjP[s]) break;
9476185916SToby Isaac           }
9576185916SToby Isaac           if (s < numAdjP) continue;
969566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(section, qAdj, &qAdjDof));
979566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetConstraintDof(section, qAdj, &qAdjCDof));
9876185916SToby Isaac           iNew += qAdjDof - qAdjCDof;
9976185916SToby Isaac         }
1009566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(adjSec, p, iNew));
10176185916SToby Isaac       }
10276185916SToby Isaac     }
10376185916SToby Isaac 
1049566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(adjSec));
1059566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(adjSec, &adjSize));
1069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(adjSize, &adj));
10776185916SToby Isaac 
10876185916SToby Isaac     for (p = pStart; p < pEnd; p++) {
10976185916SToby Isaac       PetscInt iDof, iOff, i, r, s, aOff, aOffOrig, aDof, numAdjP = PETSC_DETERMINE;
11076185916SToby Isaac 
1119566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(inverseSec, p, &iDof));
11276185916SToby Isaac       if (!iDof) continue;
1139566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(inverseSec, p, &iOff));
1149566063dSJacob Faibussowitsch       PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, PETSC_TRUE, &numAdjP, &tmpAdjP));
1159566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(adjSec, p, &aDof));
1169566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(adjSec, p, &aOff));
11776185916SToby Isaac       aOffOrig = aOff;
11876185916SToby Isaac       for (i = 0; i < iDof; i++) {
11976185916SToby Isaac         PetscInt qAdj, qAdjDof, qAdjCDof, qAdjOff, nd, numAdjQ = PETSC_DETERMINE;
12076185916SToby Isaac 
12176185916SToby Isaac         q = inverse[iOff + i];
1229566063dSJacob Faibussowitsch         PetscCall(DMPlexGetAdjacency_Internal(dm, q, useCone, useClosure, PETSC_TRUE, &numAdjQ, &tmpAdjQ));
12376185916SToby Isaac         for (r = 0; r < numAdjQ; r++) {
12476185916SToby Isaac           qAdj = tmpAdjQ[r];
12576185916SToby Isaac           if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
12676185916SToby Isaac           for (s = 0; s < numAdjP; s++) {
12776185916SToby Isaac             if (qAdj == tmpAdjP[s]) break;
12876185916SToby Isaac           }
12976185916SToby Isaac           if (s < numAdjP) continue;
1309566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(section, qAdj, &qAdjDof));
1319566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetConstraintDof(section, qAdj, &qAdjCDof));
1329566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetOffset(sectionGlobal, qAdj, &qAdjOff));
133ad540459SPierre Jolivet           for (nd = 0; nd < qAdjDof - qAdjCDof; ++nd) adj[aOff++] = (qAdjOff < 0 ? -(qAdjOff + 1) : qAdjOff) + nd;
13476185916SToby Isaac         }
13576185916SToby Isaac       }
1369566063dSJacob Faibussowitsch       PetscCall(PetscSortRemoveDupsInt(&aDof, &adj[aOffOrig]));
1379566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(adjSec, p, aDof));
13876185916SToby Isaac     }
13976185916SToby Isaac     *anchorAdj = adj;
14076185916SToby Isaac 
14176185916SToby Isaac     /* clean up */
1429566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&inverseSec));
1439566063dSJacob Faibussowitsch     PetscCall(PetscFree(inverse));
1449566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmpAdjP));
1459566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmpAdjQ));
1469371c9d4SSatish Balay   } else {
14776185916SToby Isaac     *anchorAdj = NULL;
1489566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(adjSec));
14976185916SToby Isaac   }
15076185916SToby Isaac   *anchorSectionAdj = adjSec;
15176185916SToby Isaac   PetscFunctionReturn(0);
15276185916SToby Isaac }
15376185916SToby Isaac 
154*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateAdjacencySection_Static(DM dm, PetscInt bs, PetscSF sfDof, PetscBool useCone, PetscBool useClosure, PetscBool useAnchors, PetscSection *sA, PetscInt **colIdx)
155*d71ae5a4SJacob Faibussowitsch {
156cb1e1211SMatthew G Knepley   MPI_Comm           comm;
157a9fb6443SMatthew G. Knepley   PetscMPIInt        size;
158a9fb6443SMatthew G. Knepley   PetscBool          doCommLocal, doComm, debug = PETSC_FALSE;
159a9fb6443SMatthew G. Knepley   PetscSF            sf, sfAdj;
160e101f074SMatthew G. Knepley   PetscSection       section, sectionGlobal, leafSectionAdj, rootSectionAdj, sectionAdj, anchorSectionAdj;
161a9fb6443SMatthew G. Knepley   PetscInt           nroots, nleaves, l, p, r;
162cb1e1211SMatthew G Knepley   const PetscInt    *leaves;
163cb1e1211SMatthew G Knepley   const PetscSFNode *remotes;
164cb1e1211SMatthew G Knepley   PetscInt           dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols;
1658821704fSMatthew G. Knepley   PetscInt          *tmpAdj = NULL, *adj, *rootAdj, *anchorAdj = NULL, *cols, *remoteOffsets;
16670034214SMatthew G. Knepley   PetscInt           adjSize;
167cb1e1211SMatthew G Knepley 
168cb1e1211SMatthew G Knepley   PetscFunctionBegin;
1699566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_view_preallocation", &debug, NULL));
1719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1729566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1739566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
1749566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
1759566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
1769566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
17747a6104aSMatthew G. Knepley   doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE;
1781c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm));
179cb1e1211SMatthew G Knepley   /* Create section for dof adjacency (dof ==> # adj dof) */
1809566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1819566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &numDof));
1829566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &leafSectionAdj));
1839566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(leafSectionAdj, 0, numDof));
1849566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &rootSectionAdj));
1859566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(rootSectionAdj, 0, numDof));
186cb1e1211SMatthew G Knepley   /*   Fill in the ghost dofs on the interface */
1879566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes));
188cb1e1211SMatthew G Knepley   /*
189964bf7afSMatthew G. Knepley    section        - maps points to (# dofs, local dofs)
190964bf7afSMatthew G. Knepley    sectionGlobal  - maps points to (# dofs, global dofs)
191964bf7afSMatthew G. Knepley    leafSectionAdj - maps unowned local dofs to # adj dofs
192964bf7afSMatthew G. Knepley    rootSectionAdj - maps   owned local dofs to # adj dofs
193964bf7afSMatthew G. Knepley    adj            - adj global dofs indexed by leafSectionAdj
194964bf7afSMatthew G. Knepley    rootAdj        - adj global dofs indexed by rootSectionAdj
195964bf7afSMatthew G. Knepley    sf    - describes shared points across procs
196964bf7afSMatthew G. Knepley    sfDof - describes shared dofs across procs
197964bf7afSMatthew G. Knepley    sfAdj - describes shared adjacent dofs across procs
198cb1e1211SMatthew G Knepley    ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point.
19976185916SToby Isaac   (0). If there are point-to-point constraints, add the adjacencies of constrained points to anchors in anchorAdj
20076185916SToby Isaac        (This is done in DMPlexComputeAnchorAdjacencies())
201cb1e1211SMatthew G Knepley     1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj
202cb1e1211SMatthew G Knepley        Reduce those counts to rootSectionAdj (now redundantly counting some interface points)
203cb1e1211SMatthew G Knepley     2. Visit owned points on interface, count adjacencies placing in rootSectionAdj
204cb1e1211SMatthew G Knepley        Create sfAdj connecting rootSectionAdj and leafSectionAdj
205cb1e1211SMatthew G Knepley     3. Visit unowned points on interface, write adjacencies to adj
206cb1e1211SMatthew G Knepley        Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies)
207cb1e1211SMatthew G Knepley     4. Visit owned points on interface, write adjacencies to rootAdj
208cb1e1211SMatthew G Knepley        Remove redundancy in rootAdj
209cb1e1211SMatthew G Knepley    ** The last two traversals use transitive closure
210cb1e1211SMatthew G Knepley     5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj)
211cb1e1211SMatthew G Knepley        Allocate memory addressed by sectionAdj (cols)
212cb1e1211SMatthew G Knepley     6. Visit all owned points in the subdomain, insert dof adjacencies into cols
213cb1e1211SMatthew G Knepley    ** Knowing all the column adjacencies, check ownership and sum into dnz and onz
214cb1e1211SMatthew G Knepley   */
2159566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeAnchorAdjacencies(dm, useCone, useClosure, &anchorSectionAdj, &anchorAdj));
216cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
21776185916SToby Isaac     PetscInt dof, off, d, q, anDof;
21870034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
219cb1e1211SMatthew G Knepley 
220fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
2219566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
2229566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
2239566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
224cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
225cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
226cb1e1211SMatthew G Knepley       PetscInt       ndof, ncdof;
227cb1e1211SMatthew G Knepley 
228cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
2299566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
2309566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
23148a46eb9SPierre Jolivet       for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(leafSectionAdj, d, ndof - ncdof));
232cb1e1211SMatthew G Knepley     }
2339566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
23476185916SToby Isaac     if (anDof) {
23548a46eb9SPierre Jolivet       for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(leafSectionAdj, d, anDof));
23676185916SToby Isaac     }
237cb1e1211SMatthew G Knepley   }
2389566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSectionAdj));
239cb1e1211SMatthew G Knepley   if (debug) {
2409566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n"));
2419566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(leafSectionAdj, NULL));
242cb1e1211SMatthew G Knepley   }
243cb1e1211SMatthew G Knepley   /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */
24447a6104aSMatthew G. Knepley   if (doComm) {
2459566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
2469566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
24769c11d05SVaclav Hapla     PetscCall(PetscSectionInvalidateMaxDof_Internal(rootSectionAdj));
248cb1e1211SMatthew G Knepley   }
249cb1e1211SMatthew G Knepley   if (debug) {
2509566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n"));
2519566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
252cb1e1211SMatthew G Knepley   }
253cb1e1211SMatthew G Knepley   /* Add in local adjacency sizes for owned dofs on interface (roots) */
254cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
25576185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof;
256cb1e1211SMatthew G Knepley 
2579566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
2589566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
259cb1e1211SMatthew G Knepley     if (!dof) continue;
2609566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
261cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
2629566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
263cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
264cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
265cb1e1211SMatthew G Knepley       PetscInt       ndof, ncdof;
266cb1e1211SMatthew G Knepley 
267cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
2689566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
2699566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
27048a46eb9SPierre Jolivet       for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(rootSectionAdj, d, ndof - ncdof));
271cb1e1211SMatthew G Knepley     }
2729566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
27376185916SToby Isaac     if (anDof) {
27448a46eb9SPierre Jolivet       for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(rootSectionAdj, d, anDof));
27576185916SToby Isaac     }
276cb1e1211SMatthew G Knepley   }
2779566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(rootSectionAdj));
278cb1e1211SMatthew G Knepley   if (debug) {
2799566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n"));
2809566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
281cb1e1211SMatthew G Knepley   }
282cb1e1211SMatthew G Knepley   /* Create adj SF based on dof SF */
2839566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets));
2849566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj));
2859566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
2866cded2eaSMatthew G. Knepley   if (debug && size > 1) {
2879566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency SF for Preallocation:\n"));
2889566063dSJacob Faibussowitsch     PetscCall(PetscSFView(sfAdj, NULL));
289cb1e1211SMatthew G Knepley   }
290cb1e1211SMatthew G Knepley   /* Create leaf adjacency */
2919566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSectionAdj));
2929566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(leafSectionAdj, &adjSize));
2939566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(adjSize, &adj));
294cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
29576185916SToby Isaac     PetscInt dof, off, d, q, anDof, anOff;
29670034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
297cb1e1211SMatthew G Knepley 
298fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
2999566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
3009566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
3019566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
3029566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
3039566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
304cb1e1211SMatthew G Knepley     for (d = off; d < off + dof; ++d) {
305cb1e1211SMatthew G Knepley       PetscInt aoff, i = 0;
306cb1e1211SMatthew G Knepley 
3079566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(leafSectionAdj, d, &aoff));
308cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
309cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
310cb1e1211SMatthew G Knepley         PetscInt       ndof, ncdof, ngoff, nd;
311cb1e1211SMatthew G Knepley 
312cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
3139566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
3149566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
3159566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
316cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof - ncdof; ++nd) {
317cb1e1211SMatthew G Knepley           adj[aoff + i] = (ngoff < 0 ? -(ngoff + 1) : ngoff) + nd;
318cb1e1211SMatthew G Knepley           ++i;
319cb1e1211SMatthew G Knepley         }
320cb1e1211SMatthew G Knepley       }
32176185916SToby Isaac       for (q = 0; q < anDof; q++) {
32276185916SToby Isaac         adj[aoff + i] = anchorAdj[anOff + q];
32376185916SToby Isaac         ++i;
32476185916SToby Isaac       }
325cb1e1211SMatthew G Knepley     }
326cb1e1211SMatthew G Knepley   }
327cb1e1211SMatthew G Knepley   /* Debugging */
328cb1e1211SMatthew G Knepley   if (debug) {
329cb1e1211SMatthew G Knepley     IS tmp;
3309566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Leaf adjacency indices\n"));
3319566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp));
3329566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
3339566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
334cb1e1211SMatthew G Knepley   }
335543482b8SMatthew G. Knepley   /* Gather adjacent indices to root */
3369566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(rootSectionAdj, &adjSize));
3379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(adjSize, &rootAdj));
338cb1e1211SMatthew G Knepley   for (r = 0; r < adjSize; ++r) rootAdj[r] = -1;
33947a6104aSMatthew G. Knepley   if (doComm) {
340543482b8SMatthew G. Knepley     const PetscInt *indegree;
341543482b8SMatthew G. Knepley     PetscInt       *remoteadj, radjsize = 0;
342543482b8SMatthew G. Knepley 
3439566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sfAdj, &indegree));
3449566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sfAdj, &indegree));
345543482b8SMatthew G. Knepley     for (p = 0; p < adjSize; ++p) radjsize += indegree[p];
3469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(radjsize, &remoteadj));
3479566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj));
3489566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherEnd(sfAdj, MPIU_INT, adj, remoteadj));
3496dba2905SMatthew G. Knepley     for (p = 0, l = 0, r = 0; p < adjSize; ++p, l = PetscMax(p, l + indegree[p - 1])) {
350543482b8SMatthew G. Knepley       PetscInt s;
3516dba2905SMatthew G. Knepley       for (s = 0; s < indegree[p]; ++s, ++r) rootAdj[l + s] = remoteadj[r];
352543482b8SMatthew G. Knepley     }
35363a3b9bcSJacob Faibussowitsch     PetscCheck(r == radjsize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %" PetscInt_FMT " != %" PetscInt_FMT, r, radjsize);
35463a3b9bcSJacob Faibussowitsch     PetscCheck(l == adjSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %" PetscInt_FMT " != %" PetscInt_FMT, l, adjSize);
3559566063dSJacob Faibussowitsch     PetscCall(PetscFree(remoteadj));
356cb1e1211SMatthew G Knepley   }
3579566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfAdj));
3589566063dSJacob Faibussowitsch   PetscCall(PetscFree(adj));
359cb1e1211SMatthew G Knepley   /* Debugging */
360cb1e1211SMatthew G Knepley   if (debug) {
361cb1e1211SMatthew G Knepley     IS tmp;
3629566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices after gather\n"));
3639566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
3649566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
3659566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
366cb1e1211SMatthew G Knepley   }
367cb1e1211SMatthew G Knepley   /* Add in local adjacency indices for owned dofs on interface (roots) */
368cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
36976185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff;
370cb1e1211SMatthew G Knepley 
3719566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
3729566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
373cb1e1211SMatthew G Knepley     if (!dof) continue;
3749566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
375cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
3769566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
3779566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
3789566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
379cb1e1211SMatthew G Knepley     for (d = off; d < off + dof; ++d) {
380cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i;
381cb1e1211SMatthew G Knepley 
3829566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof));
3839566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
384cb1e1211SMatthew G Knepley       i = adof - 1;
38576185916SToby Isaac       for (q = 0; q < anDof; q++) {
38676185916SToby Isaac         rootAdj[aoff + i] = anchorAdj[anOff + q];
38776185916SToby Isaac         --i;
38876185916SToby Isaac       }
389cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
390cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
391cb1e1211SMatthew G Knepley         PetscInt       ndof, ncdof, ngoff, nd;
392cb1e1211SMatthew G Knepley 
393cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
3949566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
3959566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
3969566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
397cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof - ncdof; ++nd) {
398cb1e1211SMatthew G Knepley           rootAdj[aoff + i] = ngoff < 0 ? -(ngoff + 1) + nd : ngoff + nd;
399cb1e1211SMatthew G Knepley           --i;
400cb1e1211SMatthew G Knepley         }
401cb1e1211SMatthew G Knepley       }
402cb1e1211SMatthew G Knepley     }
403cb1e1211SMatthew G Knepley   }
404cb1e1211SMatthew G Knepley   /* Debugging */
405cb1e1211SMatthew G Knepley   if (debug) {
406cb1e1211SMatthew G Knepley     IS tmp;
4079566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices\n"));
4089566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
4099566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
4109566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
411cb1e1211SMatthew G Knepley   }
412cb1e1211SMatthew G Knepley   /* Compress indices */
4139566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(rootSectionAdj));
414cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
415cb1e1211SMatthew G Knepley     PetscInt dof, cdof, off, d;
416cb1e1211SMatthew G Knepley     PetscInt adof, aoff;
417cb1e1211SMatthew G Knepley 
4189566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4199566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
4209566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
421cb1e1211SMatthew G Knepley     if (!dof) continue;
4229566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
423cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
424cb1e1211SMatthew G Knepley     for (d = off; d < off + dof - cdof; ++d) {
4259566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof));
4269566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
4279566063dSJacob Faibussowitsch       PetscCall(PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]));
4289566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(rootSectionAdj, d, adof));
429cb1e1211SMatthew G Knepley     }
430cb1e1211SMatthew G Knepley   }
431cb1e1211SMatthew G Knepley   /* Debugging */
432cb1e1211SMatthew G Knepley   if (debug) {
433cb1e1211SMatthew G Knepley     IS tmp;
4349566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n"));
4359566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
4369566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices after compression\n"));
4379566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
4389566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
4399566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
440cb1e1211SMatthew G Knepley   }
441cb1e1211SMatthew G Knepley   /* Build adjacency section: Maps global indices to sets of adjacent global indices */
4429566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd));
4439566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &sectionAdj));
4449566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd));
445cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
44676185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof;
447cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
448cb1e1211SMatthew G Knepley 
4499566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4509566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
4519566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
4529566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
453cb1e1211SMatthew G Knepley     for (d = 0; d < dof - cdof; ++d) {
454cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
455cb1e1211SMatthew G Knepley 
4569566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSectionAdj, off + d, &ldof));
4579566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, off + d, &rdof));
458cb1e1211SMatthew G Knepley       if (ldof > 0) {
459cb1e1211SMatthew G Knepley         /* We do not own this point */
460cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
4619566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetDof(sectionAdj, goff + d, rdof));
462cb1e1211SMatthew G Knepley       } else {
463cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
464cb1e1211SMatthew G Knepley       }
465cb1e1211SMatthew G Knepley     }
466cb1e1211SMatthew G Knepley     if (found) continue;
4679566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4689566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
4699566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
470cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
471cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
472cb1e1211SMatthew G Knepley       PetscInt       ndof, ncdof, noff;
473cb1e1211SMatthew G Knepley 
474cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
4759566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
4769566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
4779566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, padj, &noff));
47848a46eb9SPierre Jolivet       for (d = goff; d < goff + dof - cdof; ++d) PetscCall(PetscSectionAddDof(sectionAdj, d, ndof - ncdof));
479cb1e1211SMatthew G Knepley     }
4809566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
48176185916SToby Isaac     if (anDof) {
48248a46eb9SPierre Jolivet       for (d = goff; d < goff + dof - cdof; ++d) PetscCall(PetscSectionAddDof(sectionAdj, d, anDof));
48376185916SToby Isaac     }
484cb1e1211SMatthew G Knepley   }
4859566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(sectionAdj));
486cb1e1211SMatthew G Knepley   if (debug) {
4879566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation:\n"));
4889566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(sectionAdj, NULL));
489cb1e1211SMatthew G Knepley   }
490cb1e1211SMatthew G Knepley   /* Get adjacent indices */
4919566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(sectionAdj, &numCols));
4929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numCols, &cols));
493cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
49476185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff;
495cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
496cb1e1211SMatthew G Knepley 
4979566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4989566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
4999566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
5009566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
501cb1e1211SMatthew G Knepley     for (d = 0; d < dof - cdof; ++d) {
502cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
503cb1e1211SMatthew G Knepley 
5049566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSectionAdj, off + d, &ldof));
5059566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, off + d, &rdof));
506cb1e1211SMatthew G Knepley       if (ldof > 0) {
507cb1e1211SMatthew G Knepley         /* We do not own this point */
508cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
509cb1e1211SMatthew G Knepley         PetscInt aoff, roff;
510cb1e1211SMatthew G Knepley 
5119566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, goff + d, &aoff));
5129566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(rootSectionAdj, off + d, &roff));
5139566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&cols[aoff], &rootAdj[roff], rdof));
514cb1e1211SMatthew G Knepley       } else {
515cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
516cb1e1211SMatthew G Knepley       }
517cb1e1211SMatthew G Knepley     }
518cb1e1211SMatthew G Knepley     if (found) continue;
5199566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
5209566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
5219566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
522cb1e1211SMatthew G Knepley     for (d = goff; d < goff + dof - cdof; ++d) {
523cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i = 0;
524cb1e1211SMatthew G Knepley 
5259566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, d, &adof));
5269566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(sectionAdj, d, &aoff));
527cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
528cb1e1211SMatthew G Knepley         const PetscInt  padj = tmpAdj[q];
529cb1e1211SMatthew G Knepley         PetscInt        ndof, ncdof, ngoff, nd;
530cb1e1211SMatthew G Knepley         const PetscInt *ncind;
531cb1e1211SMatthew G Knepley 
532cb1e1211SMatthew G Knepley         /* Adjacent points may not be in the section chart */
533cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
5349566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
5359566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
5369566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintIndices(section, padj, &ncind));
5379566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
538ad540459SPierre Jolivet         for (nd = 0; nd < ndof - ncdof; ++nd, ++i) cols[aoff + i] = ngoff < 0 ? -(ngoff + 1) + nd : ngoff + nd;
539cb1e1211SMatthew G Knepley       }
540ad540459SPierre Jolivet       for (q = 0; q < anDof; q++, i++) cols[aoff + i] = anchorAdj[anOff + q];
54163a3b9bcSJacob 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);
542cb1e1211SMatthew G Knepley     }
543cb1e1211SMatthew G Knepley   }
5449566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&anchorSectionAdj));
5459566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&leafSectionAdj));
5469566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&rootSectionAdj));
5479566063dSJacob Faibussowitsch   PetscCall(PetscFree(anchorAdj));
5489566063dSJacob Faibussowitsch   PetscCall(PetscFree(rootAdj));
5499566063dSJacob Faibussowitsch   PetscCall(PetscFree(tmpAdj));
550cb1e1211SMatthew G Knepley   /* Debugging */
551cb1e1211SMatthew G Knepley   if (debug) {
552cb1e1211SMatthew G Knepley     IS tmp;
5539566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Column indices\n"));
5549566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp));
5559566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
5569566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
557cb1e1211SMatthew G Knepley   }
558a9fb6443SMatthew G. Knepley 
559a9fb6443SMatthew G. Knepley   *sA     = sectionAdj;
560a9fb6443SMatthew G. Knepley   *colIdx = cols;
561a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
562a9fb6443SMatthew G. Knepley }
563a9fb6443SMatthew G. Knepley 
564*d71ae5a4SJacob Faibussowitsch 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[])
565*d71ae5a4SJacob Faibussowitsch {
566e101f074SMatthew G. Knepley   PetscSection section;
567a9fb6443SMatthew G. Knepley   PetscInt     rStart, rEnd, r, pStart, pEnd, p;
568a9fb6443SMatthew G. Knepley 
569a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
570a9fb6443SMatthew G. Knepley   /* This loop needs to change to a loop over points, then field dofs, which means we need to look both sections */
5719566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
5721dca8a05SBarry 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);
573a9fb6443SMatthew G. Knepley   if (f >= 0 && bs == 1) {
5749566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
5759566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
576a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
577294b7009SMatthew G. Knepley       PetscInt rS, rE;
578a9fb6443SMatthew G. Knepley 
5799566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE));
580a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
581a9fb6443SMatthew G. Knepley         PetscInt numCols, cStart, c;
582a9fb6443SMatthew G. Knepley 
5839566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
5849566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
585a9fb6443SMatthew G. Knepley         for (c = cStart; c < cStart + numCols; ++c) {
586a9fb6443SMatthew G. Knepley           if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
587a9fb6443SMatthew G. Knepley             ++dnz[r - rStart];
588a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++dnzu[r - rStart];
589a9fb6443SMatthew G. Knepley           } else {
590a9fb6443SMatthew G. Knepley             ++onz[r - rStart];
591a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++onzu[r - rStart];
592a9fb6443SMatthew G. Knepley           }
593a9fb6443SMatthew G. Knepley         }
594a9fb6443SMatthew G. Knepley       }
595a9fb6443SMatthew G. Knepley     }
596a9fb6443SMatthew G. Knepley   } else {
597cb1e1211SMatthew G Knepley     /* Only loop over blocks of rows */
598cb1e1211SMatthew G Knepley     for (r = rStart / bs; r < rEnd / bs; ++r) {
599cb1e1211SMatthew G Knepley       const PetscInt row = r * bs;
600cb1e1211SMatthew G Knepley       PetscInt       numCols, cStart, c;
601cb1e1211SMatthew G Knepley 
6029566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, row, &numCols));
6039566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(sectionAdj, row, &cStart));
604cb1e1211SMatthew G Knepley       for (c = cStart; c < cStart + numCols; ++c) {
605e7bcfa36SMatthew G. Knepley         if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
606e7bcfa36SMatthew G. Knepley           ++dnz[r - rStart / bs];
607e7bcfa36SMatthew G. Knepley           if (cols[c] >= row) ++dnzu[r - rStart / bs];
608cb1e1211SMatthew G Knepley         } else {
609e7bcfa36SMatthew G. Knepley           ++onz[r - rStart / bs];
610e7bcfa36SMatthew G. Knepley           if (cols[c] >= row) ++onzu[r - rStart / bs];
611cb1e1211SMatthew G Knepley         }
612cb1e1211SMatthew G Knepley       }
613cb1e1211SMatthew G Knepley     }
614a9fb6443SMatthew G. Knepley     for (r = 0; r < (rEnd - rStart) / bs; ++r) {
615cb1e1211SMatthew G Knepley       dnz[r] /= bs;
616cb1e1211SMatthew G Knepley       onz[r] /= bs;
617cb1e1211SMatthew G Knepley       dnzu[r] /= bs;
618cb1e1211SMatthew G Knepley       onzu[r] /= bs;
619cb1e1211SMatthew G Knepley     }
620cb1e1211SMatthew G Knepley   }
621a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
622a9fb6443SMatthew G. Knepley }
623a9fb6443SMatthew G. Knepley 
624*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexFillMatrix_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], Mat A)
625*d71ae5a4SJacob Faibussowitsch {
626e101f074SMatthew G. Knepley   PetscSection section;
627a9fb6443SMatthew G. Knepley   PetscScalar *values;
628a9fb6443SMatthew G. Knepley   PetscInt     rStart, rEnd, r, pStart, pEnd, p, len, maxRowLen = 0;
629a9fb6443SMatthew G. Knepley 
630a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
6319566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
632a9fb6443SMatthew G. Knepley   for (r = rStart; r < rEnd; ++r) {
6339566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(sectionAdj, r, &len));
634a9fb6443SMatthew G. Knepley     maxRowLen = PetscMax(maxRowLen, len);
635a9fb6443SMatthew G. Knepley   }
6369566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(maxRowLen, &values));
637a9fb6443SMatthew G. Knepley   if (f >= 0 && bs == 1) {
6389566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
6399566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
640a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
641294b7009SMatthew G. Knepley       PetscInt rS, rE;
642a9fb6443SMatthew G. Knepley 
6439566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE));
644a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
6457e01ee02SMatthew G. Knepley         PetscInt numCols, cStart;
646a9fb6443SMatthew G. Knepley 
6479566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
6489566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
6499566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
650a9fb6443SMatthew G. Knepley       }
651a9fb6443SMatthew G. Knepley     }
652a9fb6443SMatthew G. Knepley   } else {
653a9fb6443SMatthew G. Knepley     for (r = rStart; r < rEnd; ++r) {
654a9fb6443SMatthew G. Knepley       PetscInt numCols, cStart;
655a9fb6443SMatthew G. Knepley 
6569566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
6579566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
6589566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
659a9fb6443SMatthew G. Knepley     }
660a9fb6443SMatthew G. Knepley   }
6619566063dSJacob Faibussowitsch   PetscCall(PetscFree(values));
662a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
663a9fb6443SMatthew G. Knepley }
664a9fb6443SMatthew G. Knepley 
66525e402d2SMatthew G. Knepley /*@C
66625e402d2SMatthew G. Knepley   DMPlexPreallocateOperator - Calculate the matrix nonzero pattern based upon the information in the DM,
66725e402d2SMatthew G. Knepley   the PetscDS it contains, and the default PetscSection.
66825e402d2SMatthew G. Knepley 
66925e402d2SMatthew G. Knepley   Collective
67025e402d2SMatthew G. Knepley 
6714165533cSJose E. Roman   Input Parameters:
67225e402d2SMatthew G. Knepley + dm   - The DMPlex
67325e402d2SMatthew G. Knepley . bs   - The matrix blocksize
67425e402d2SMatthew G. Knepley . dnz  - An array to hold the number of nonzeros in the diagonal block
67525e402d2SMatthew G. Knepley . onz  - An array to hold the number of nonzeros in the off-diagonal block
67625e402d2SMatthew G. Knepley . dnzu - An array to hold the number of nonzeros in the upper triangle of the diagonal block
67725e402d2SMatthew G. Knepley . onzu - An array to hold the number of nonzeros in the upper triangle of the off-diagonal block
678a2b725a8SWilliam Gropp - fillMatrix - If PETSC_TRUE, fill the matrix with zeros
67925e402d2SMatthew G. Knepley 
6804165533cSJose E. Roman   Output Parameter:
68125e402d2SMatthew G. Knepley . A - The preallocated matrix
68225e402d2SMatthew G. Knepley 
68325e402d2SMatthew G. Knepley   Level: advanced
68425e402d2SMatthew G. Knepley 
685db781477SPatrick Sanan .seealso: `DMCreateMatrix()`
68625e402d2SMatthew G. Knepley @*/
687*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix)
688*d71ae5a4SJacob Faibussowitsch {
689a9fb6443SMatthew G. Knepley   MPI_Comm     comm;
690a9fb6443SMatthew G. Knepley   PetscDS      prob;
691a9fb6443SMatthew G. Knepley   MatType      mtype;
692a9fb6443SMatthew G. Knepley   PetscSF      sf, sfDof;
693e101f074SMatthew G. Knepley   PetscSection section;
694a9fb6443SMatthew G. Knepley   PetscInt    *remoteOffsets;
695a9fb6443SMatthew G. Knepley   PetscSection sectionAdj[4] = {NULL, NULL, NULL, NULL};
696a9fb6443SMatthew G. Knepley   PetscInt    *cols[4]       = {NULL, NULL, NULL, NULL};
697a9fb6443SMatthew G. Knepley   PetscBool    useCone, useClosure;
6987e01ee02SMatthew G. Knepley   PetscInt     Nf, f, idx, locRows;
699a9fb6443SMatthew G. Knepley   PetscLayout  rLayout;
700a9fb6443SMatthew G. Knepley   PetscBool    isSymBlock, isSymSeqBlock, isSymMPIBlock, debug = PETSC_FALSE;
7016cded2eaSMatthew G. Knepley   PetscMPIInt  size;
702a9fb6443SMatthew G. Knepley 
703a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
704a9fb6443SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
705064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(A, MAT_CLASSID, 7);
706dadcf809SJacob Faibussowitsch   if (dnz) PetscValidIntPointer(dnz, 3);
707dadcf809SJacob Faibussowitsch   if (onz) PetscValidIntPointer(onz, 4);
708dadcf809SJacob Faibussowitsch   if (dnzu) PetscValidIntPointer(dnzu, 5);
709dadcf809SJacob Faibussowitsch   if (onzu) PetscValidIntPointer(onzu, 6);
7109566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
7119566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
7129566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
7139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_view_preallocation", &debug, NULL));
7149566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
7159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
7169566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_Preallocate, dm, 0, 0, 0));
717a9fb6443SMatthew G. Knepley   /* Create dof SF based on point SF */
718a9fb6443SMatthew G. Knepley   if (debug) {
719e101f074SMatthew G. Knepley     PetscSection section, sectionGlobal;
720a9fb6443SMatthew G. Knepley     PetscSF      sf;
721a9fb6443SMatthew G. Knepley 
7229566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
7239566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
7249566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
7259566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Input Section for Preallocation:\n"));
7269566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(section, NULL));
7279566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Input Global Section for Preallocation:\n"));
7289566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(sectionGlobal, NULL));
7296cded2eaSMatthew G. Knepley     if (size > 1) {
7309566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Input SF for Preallocation:\n"));
7319566063dSJacob Faibussowitsch       PetscCall(PetscSFView(sf, NULL));
732a9fb6443SMatthew G. Knepley     }
7336cded2eaSMatthew G. Knepley   }
7349566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets));
7359566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof));
7369566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
7376cded2eaSMatthew G. Knepley   if (debug && size > 1) {
7389566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Dof SF for Preallocation:\n"));
7399566063dSJacob Faibussowitsch     PetscCall(PetscSFView(sfDof, NULL));
740a9fb6443SMatthew G. Knepley   }
741a9fb6443SMatthew G. Knepley   /* Create allocation vectors from adjacency graph */
7429566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &locRows, NULL));
7439566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
7449566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
7459566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
7469566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
747a9fb6443SMatthew G. Knepley   /* There are 4 types of adjacency */
7489566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
749a9fb6443SMatthew G. Knepley   if (Nf < 1 || bs > 1) {
7509566063dSJacob Faibussowitsch     PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
751a9fb6443SMatthew G. Knepley     idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7529566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]));
7539566063dSJacob Faibussowitsch     PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
754a9fb6443SMatthew G. Knepley   } else {
755a9fb6443SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
7569566063dSJacob Faibussowitsch       PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
757a9fb6443SMatthew G. Knepley       idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7589566063dSJacob Faibussowitsch       if (!sectionAdj[idx]) PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]));
7599566063dSJacob Faibussowitsch       PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
760a9fb6443SMatthew G. Knepley     }
761a9fb6443SMatthew G. Knepley   }
7629566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfDof));
763cb1e1211SMatthew G Knepley   /* Set matrix pattern */
7649566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu));
7659566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
76689545effSMatthew G. Knepley   /* Check for symmetric storage */
7679566063dSJacob Faibussowitsch   PetscCall(MatGetType(A, &mtype));
7689566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATSBAIJ, &isSymBlock));
7699566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock));
7709566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock));
7719566063dSJacob Faibussowitsch   if (isSymBlock || isSymSeqBlock || isSymMPIBlock) PetscCall(MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
772cb1e1211SMatthew G Knepley   /* Fill matrix with zeros */
773cb1e1211SMatthew G Knepley   if (fillMatrix) {
774a9fb6443SMatthew G. Knepley     if (Nf < 1 || bs > 1) {
7759566063dSJacob Faibussowitsch       PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
776a9fb6443SMatthew G. Knepley       idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7779566063dSJacob Faibussowitsch       PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], A));
778a9fb6443SMatthew G. Knepley     } else {
779a9fb6443SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
7809566063dSJacob Faibussowitsch         PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
781a9fb6443SMatthew G. Knepley         idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7829566063dSJacob Faibussowitsch         PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], A));
783cb1e1211SMatthew G Knepley       }
784cb1e1211SMatthew G Knepley     }
7859566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
7869566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
787cb1e1211SMatthew G Knepley   }
7889566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
7899371c9d4SSatish Balay   for (idx = 0; idx < 4; ++idx) {
7909371c9d4SSatish Balay     PetscCall(PetscSectionDestroy(&sectionAdj[idx]));
7919371c9d4SSatish Balay     PetscCall(PetscFree(cols[idx]));
7929371c9d4SSatish Balay   }
7939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_Preallocate, dm, 0, 0, 0));
794cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
795cb1e1211SMatthew G Knepley }
796cb1e1211SMatthew G Knepley 
797cb1e1211SMatthew G Knepley #if 0
798cb1e1211SMatthew 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)
799cb1e1211SMatthew G Knepley {
800cb1e1211SMatthew G Knepley   PetscInt       *tmpClosure,*tmpAdj,*visits;
801cb1e1211SMatthew G Knepley   PetscInt        c,cStart,cEnd,pStart,pEnd;
802cb1e1211SMatthew G Knepley 
803cb1e1211SMatthew G Knepley   PetscFunctionBegin;
8049566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
8059566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
8069566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
807cb1e1211SMatthew G Knepley 
808e7b80042SMatthew G. Knepley   maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1));
809cb1e1211SMatthew G Knepley 
8109566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
811cb1e1211SMatthew G Knepley   npoints = pEnd - pStart;
812cb1e1211SMatthew G Knepley 
8139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits));
8149566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(lvisits,pEnd-pStart));
8159566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(visits,pEnd-pStart));
8169566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
817cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
818cb1e1211SMatthew G Knepley     PetscInt *support = tmpClosure;
8199566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support));
820cb1e1211SMatthew G Knepley     for (p=0; p<supportSize; p++) lvisits[support[p]]++;
821cb1e1211SMatthew G Knepley   }
8229566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM));
8239566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd  (sf,MPIU_INT,lvisits,visits,MPI_SUM));
8249566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits,MPI_REPLACE));
8259566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd  (sf,MPIU_INT,visits,lvisits));
826cb1e1211SMatthew G Knepley 
8279566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks());
828cb1e1211SMatthew G Knepley 
8299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner));
830cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
8319566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(cellmat,maxClosureSize*maxClosureSize));
832cb1e1211SMatthew G Knepley     /*
833cb1e1211SMatthew G Knepley      Depth-first walk of transitive closure.
834cb1e1211SMatthew 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.
835cb1e1211SMatthew G Knepley      This contribution is added to dnz if owning ranks of p and q match, to onz otherwise.
836cb1e1211SMatthew G Knepley      */
837cb1e1211SMatthew G Knepley   }
838cb1e1211SMatthew G Knepley 
8399566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM));
8409566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd  (sf,MPIU_INT,lonz,onz,MPI_SUM));
841cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
842cb1e1211SMatthew G Knepley }
843cb1e1211SMatthew G Knepley #endif
844