xref: /petsc/src/dm/impls/plex/plexpreallocate.c (revision 1c2dc1cbabe5212f80cf309c4ca5a26f0cadc660)
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() */
7e101f074SMatthew G. Knepley static PetscErrorCode DMPlexComputeAnchorAdjacencies(DM dm, PetscBool useCone, PetscBool useClosure, PetscSection *anchorSectionAdj, PetscInt *anchorAdj[])
876185916SToby Isaac {
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));
13376185916SToby Isaac           for (nd = 0; nd < qAdjDof-qAdjCDof; ++nd) {
13476185916SToby Isaac             adj[aOff++] = (qAdjOff < 0 ? -(qAdjOff+1) : qAdjOff) + nd;
13576185916SToby Isaac           }
13676185916SToby Isaac         }
13776185916SToby Isaac       }
1389566063dSJacob Faibussowitsch       PetscCall(PetscSortRemoveDupsInt(&aDof,&adj[aOffOrig]));
1399566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(adjSec,p,aDof));
14076185916SToby Isaac     }
14176185916SToby Isaac     *anchorAdj = adj;
14276185916SToby Isaac 
14376185916SToby Isaac     /* clean up */
1449566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&inverseSec));
1459566063dSJacob Faibussowitsch     PetscCall(PetscFree(inverse));
1469566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmpAdjP));
1479566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmpAdjQ));
14876185916SToby Isaac   }
14976185916SToby Isaac   else {
15076185916SToby Isaac     *anchorAdj = NULL;
1519566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(adjSec));
15276185916SToby Isaac   }
15376185916SToby Isaac   *anchorSectionAdj = adjSec;
15476185916SToby Isaac   PetscFunctionReturn(0);
15576185916SToby Isaac }
15676185916SToby Isaac 
15725e402d2SMatthew G. Knepley static PetscErrorCode DMPlexCreateAdjacencySection_Static(DM dm, PetscInt bs, PetscSF sfDof, PetscBool useCone, PetscBool useClosure, PetscBool useAnchors, PetscSection *sA, PetscInt **colIdx)
158cb1e1211SMatthew G Knepley {
159cb1e1211SMatthew G Knepley   MPI_Comm           comm;
160a9fb6443SMatthew G. Knepley   PetscMPIInt        size;
161a9fb6443SMatthew G. Knepley   PetscBool          doCommLocal, doComm, debug = PETSC_FALSE;
162a9fb6443SMatthew G. Knepley   PetscSF            sf, sfAdj;
163e101f074SMatthew G. Knepley   PetscSection       section, sectionGlobal, leafSectionAdj, rootSectionAdj, sectionAdj, anchorSectionAdj;
164a9fb6443SMatthew G. Knepley   PetscInt           nroots, nleaves, l, p, r;
165cb1e1211SMatthew G Knepley   const PetscInt    *leaves;
166cb1e1211SMatthew G Knepley   const PetscSFNode *remotes;
167cb1e1211SMatthew G Knepley   PetscInt           dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols;
1688821704fSMatthew G. Knepley   PetscInt          *tmpAdj = NULL, *adj, *rootAdj, *anchorAdj = NULL, *cols, *remoteOffsets;
16970034214SMatthew G. Knepley   PetscInt           adjSize;
170cb1e1211SMatthew G Knepley 
171cb1e1211SMatthew G Knepley   PetscFunctionBegin;
1729566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
1739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL, "-dm_view_preallocation", &debug, NULL));
1749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1759566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1769566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
1779566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
1789566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
1799566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
18047a6104aSMatthew G. Knepley   doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE;
181*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm));
182cb1e1211SMatthew G Knepley   /* Create section for dof adjacency (dof ==> # adj dof) */
1839566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1849566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &numDof));
1859566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &leafSectionAdj));
1869566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(leafSectionAdj, 0, numDof));
1879566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &rootSectionAdj));
1889566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(rootSectionAdj, 0, numDof));
189cb1e1211SMatthew G Knepley   /*   Fill in the ghost dofs on the interface */
1909566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes));
191cb1e1211SMatthew G Knepley   /*
192964bf7afSMatthew G. Knepley    section        - maps points to (# dofs, local dofs)
193964bf7afSMatthew G. Knepley    sectionGlobal  - maps points to (# dofs, global dofs)
194964bf7afSMatthew G. Knepley    leafSectionAdj - maps unowned local dofs to # adj dofs
195964bf7afSMatthew G. Knepley    rootSectionAdj - maps   owned local dofs to # adj dofs
196964bf7afSMatthew G. Knepley    adj            - adj global dofs indexed by leafSectionAdj
197964bf7afSMatthew G. Knepley    rootAdj        - adj global dofs indexed by rootSectionAdj
198964bf7afSMatthew G. Knepley    sf    - describes shared points across procs
199964bf7afSMatthew G. Knepley    sfDof - describes shared dofs across procs
200964bf7afSMatthew G. Knepley    sfAdj - describes shared adjacent dofs across procs
201cb1e1211SMatthew G Knepley    ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point.
20276185916SToby Isaac   (0). If there are point-to-point constraints, add the adjacencies of constrained points to anchors in anchorAdj
20376185916SToby Isaac        (This is done in DMPlexComputeAnchorAdjacencies())
204cb1e1211SMatthew G Knepley     1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj
205cb1e1211SMatthew G Knepley        Reduce those counts to rootSectionAdj (now redundantly counting some interface points)
206cb1e1211SMatthew G Knepley     2. Visit owned points on interface, count adjacencies placing in rootSectionAdj
207cb1e1211SMatthew G Knepley        Create sfAdj connecting rootSectionAdj and leafSectionAdj
208cb1e1211SMatthew G Knepley     3. Visit unowned points on interface, write adjacencies to adj
209cb1e1211SMatthew G Knepley        Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies)
210cb1e1211SMatthew G Knepley     4. Visit owned points on interface, write adjacencies to rootAdj
211cb1e1211SMatthew G Knepley        Remove redundancy in rootAdj
212cb1e1211SMatthew G Knepley    ** The last two traversals use transitive closure
213cb1e1211SMatthew G Knepley     5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj)
214cb1e1211SMatthew G Knepley        Allocate memory addressed by sectionAdj (cols)
215cb1e1211SMatthew G Knepley     6. Visit all owned points in the subdomain, insert dof adjacencies into cols
216cb1e1211SMatthew G Knepley    ** Knowing all the column adjacencies, check ownership and sum into dnz and onz
217cb1e1211SMatthew G Knepley   */
2189566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeAnchorAdjacencies(dm, useCone, useClosure, &anchorSectionAdj, &anchorAdj));
219cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
22076185916SToby Isaac     PetscInt dof, off, d, q, anDof;
22170034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
222cb1e1211SMatthew G Knepley 
223fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
2249566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
2259566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
2269566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
227cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
228cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
229cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
230cb1e1211SMatthew G Knepley 
231cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
2329566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
2339566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
234cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
2359566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(leafSectionAdj, d, ndof-ncdof));
236cb1e1211SMatthew G Knepley       }
237cb1e1211SMatthew G Knepley     }
2389566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
23976185916SToby Isaac     if (anDof) {
24076185916SToby Isaac       for (d = off; d < off+dof; ++d) {
2419566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(leafSectionAdj, d, anDof));
24276185916SToby Isaac       }
24376185916SToby Isaac     }
244cb1e1211SMatthew G Knepley   }
2459566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSectionAdj));
246cb1e1211SMatthew G Knepley   if (debug) {
2479566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n"));
2489566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(leafSectionAdj, NULL));
249cb1e1211SMatthew G Knepley   }
250cb1e1211SMatthew G Knepley   /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */
25147a6104aSMatthew G. Knepley   if (doComm) {
2529566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
2539566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
25469c11d05SVaclav Hapla     PetscCall(PetscSectionInvalidateMaxDof_Internal(rootSectionAdj));
255cb1e1211SMatthew G Knepley   }
256cb1e1211SMatthew G Knepley   if (debug) {
2579566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n"));
2589566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
259cb1e1211SMatthew G Knepley   }
260cb1e1211SMatthew G Knepley   /* Add in local adjacency sizes for owned dofs on interface (roots) */
261cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
26276185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof;
263cb1e1211SMatthew G Knepley 
2649566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
2659566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
266cb1e1211SMatthew G Knepley     if (!dof) continue;
2679566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
268cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
2699566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
270cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
271cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
272cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
273cb1e1211SMatthew G Knepley 
274cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
2759566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
2769566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
277cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
2789566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(rootSectionAdj, d, ndof-ncdof));
279cb1e1211SMatthew G Knepley       }
280cb1e1211SMatthew G Knepley     }
2819566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
28276185916SToby Isaac     if (anDof) {
28376185916SToby Isaac       for (d = off; d < off+dof; ++d) {
2849566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(rootSectionAdj, d, anDof));
28576185916SToby Isaac       }
28676185916SToby Isaac     }
287cb1e1211SMatthew G Knepley   }
2889566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(rootSectionAdj));
289cb1e1211SMatthew G Knepley   if (debug) {
2909566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n"));
2919566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
292cb1e1211SMatthew G Knepley   }
293cb1e1211SMatthew G Knepley   /* Create adj SF based on dof SF */
2949566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets));
2959566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj));
2969566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
2976cded2eaSMatthew G. Knepley   if (debug && size > 1) {
2989566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency SF for Preallocation:\n"));
2999566063dSJacob Faibussowitsch     PetscCall(PetscSFView(sfAdj, NULL));
300cb1e1211SMatthew G Knepley   }
301cb1e1211SMatthew G Knepley   /* Create leaf adjacency */
3029566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSectionAdj));
3039566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(leafSectionAdj, &adjSize));
3049566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(adjSize, &adj));
305cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
30676185916SToby Isaac     PetscInt dof, off, d, q, anDof, anOff;
30770034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
308cb1e1211SMatthew G Knepley 
309fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
3109566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
3119566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
3129566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
3139566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
3149566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
315cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
316cb1e1211SMatthew G Knepley       PetscInt aoff, i = 0;
317cb1e1211SMatthew G Knepley 
3189566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(leafSectionAdj, d, &aoff));
319cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
320cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
321cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
322cb1e1211SMatthew G Knepley 
323cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
3249566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
3259566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
3269566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
327cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
328cb1e1211SMatthew G Knepley           adj[aoff+i] = (ngoff < 0 ? -(ngoff+1) : ngoff) + nd;
329cb1e1211SMatthew G Knepley           ++i;
330cb1e1211SMatthew G Knepley         }
331cb1e1211SMatthew G Knepley       }
33276185916SToby Isaac       for (q = 0; q < anDof; q++) {
33376185916SToby Isaac         adj[aoff+i] = anchorAdj[anOff+q];
33476185916SToby Isaac         ++i;
33576185916SToby Isaac       }
336cb1e1211SMatthew G Knepley     }
337cb1e1211SMatthew G Knepley   }
338cb1e1211SMatthew G Knepley   /* Debugging */
339cb1e1211SMatthew G Knepley   if (debug) {
340cb1e1211SMatthew G Knepley     IS tmp;
3419566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Leaf adjacency indices\n"));
3429566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp));
3439566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
3449566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
345cb1e1211SMatthew G Knepley   }
346543482b8SMatthew G. Knepley   /* Gather adjacent indices to root */
3479566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(rootSectionAdj, &adjSize));
3489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(adjSize, &rootAdj));
349cb1e1211SMatthew G Knepley   for (r = 0; r < adjSize; ++r) rootAdj[r] = -1;
35047a6104aSMatthew G. Knepley   if (doComm) {
351543482b8SMatthew G. Knepley     const PetscInt *indegree;
352543482b8SMatthew G. Knepley     PetscInt       *remoteadj, radjsize = 0;
353543482b8SMatthew G. Knepley 
3549566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sfAdj, &indegree));
3559566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sfAdj, &indegree));
356543482b8SMatthew G. Knepley     for (p = 0; p < adjSize; ++p) radjsize += indegree[p];
3579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(radjsize, &remoteadj));
3589566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj));
3599566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherEnd(sfAdj, MPIU_INT, adj, remoteadj));
3606dba2905SMatthew G. Knepley     for (p = 0, l = 0, r = 0; p < adjSize; ++p, l = PetscMax(p, l + indegree[p-1])) {
361543482b8SMatthew G. Knepley       PetscInt s;
3626dba2905SMatthew G. Knepley       for (s = 0; s < indegree[p]; ++s, ++r) rootAdj[l+s] = remoteadj[r];
363543482b8SMatthew G. Knepley     }
3642c71b3e2SJacob Faibussowitsch     PetscCheckFalse(r != radjsize,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", r, radjsize);
3652c71b3e2SJacob Faibussowitsch     PetscCheckFalse(l != adjSize,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", l, adjSize);
3669566063dSJacob Faibussowitsch     PetscCall(PetscFree(remoteadj));
367cb1e1211SMatthew G Knepley   }
3689566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfAdj));
3699566063dSJacob Faibussowitsch   PetscCall(PetscFree(adj));
370cb1e1211SMatthew G Knepley   /* Debugging */
371cb1e1211SMatthew G Knepley   if (debug) {
372cb1e1211SMatthew G Knepley     IS tmp;
3739566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices after gather\n"));
3749566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
3759566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
3769566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
377cb1e1211SMatthew G Knepley   }
378cb1e1211SMatthew G Knepley   /* Add in local adjacency indices for owned dofs on interface (roots) */
379cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
38076185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff;
381cb1e1211SMatthew G Knepley 
3829566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
3839566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
384cb1e1211SMatthew G Knepley     if (!dof) continue;
3859566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
386cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
3879566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
3889566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
3899566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
390cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
391cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i;
392cb1e1211SMatthew G Knepley 
3939566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof));
3949566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
395cb1e1211SMatthew G Knepley       i    = adof-1;
39676185916SToby Isaac       for (q = 0; q < anDof; q++) {
39776185916SToby Isaac         rootAdj[aoff+i] = anchorAdj[anOff+q];
39876185916SToby Isaac         --i;
39976185916SToby Isaac       }
400cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
401cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
402cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
403cb1e1211SMatthew G Knepley 
404cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
4059566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
4069566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
4079566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
408cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
409cb1e1211SMatthew G Knepley           rootAdj[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
410cb1e1211SMatthew G Knepley           --i;
411cb1e1211SMatthew G Knepley         }
412cb1e1211SMatthew G Knepley       }
413cb1e1211SMatthew G Knepley     }
414cb1e1211SMatthew G Knepley   }
415cb1e1211SMatthew G Knepley   /* Debugging */
416cb1e1211SMatthew G Knepley   if (debug) {
417cb1e1211SMatthew G Knepley     IS tmp;
4189566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices\n"));
4199566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
4209566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
4219566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
422cb1e1211SMatthew G Knepley   }
423cb1e1211SMatthew G Knepley   /* Compress indices */
4249566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(rootSectionAdj));
425cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
426cb1e1211SMatthew G Knepley     PetscInt dof, cdof, off, d;
427cb1e1211SMatthew G Knepley     PetscInt adof, aoff;
428cb1e1211SMatthew G Knepley 
4299566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4309566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
4319566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
432cb1e1211SMatthew G Knepley     if (!dof) continue;
4339566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
434cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
435cb1e1211SMatthew G Knepley     for (d = off; d < off+dof-cdof; ++d) {
4369566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof));
4379566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
4389566063dSJacob Faibussowitsch       PetscCall(PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]));
4399566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(rootSectionAdj, d, adof));
440cb1e1211SMatthew G Knepley     }
441cb1e1211SMatthew G Knepley   }
442cb1e1211SMatthew G Knepley   /* Debugging */
443cb1e1211SMatthew G Knepley   if (debug) {
444cb1e1211SMatthew G Knepley     IS tmp;
4459566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n"));
4469566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
4479566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices after compression\n"));
4489566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
4499566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
4509566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
451cb1e1211SMatthew G Knepley   }
452cb1e1211SMatthew G Knepley   /* Build adjacency section: Maps global indices to sets of adjacent global indices */
4539566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd));
4549566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &sectionAdj));
4559566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd));
456cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
45776185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof;
458cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
459cb1e1211SMatthew G Knepley 
4609566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4619566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
4629566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
4639566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
464cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
465cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
466cb1e1211SMatthew G Knepley 
4679566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSectionAdj, off+d, &ldof));
4689566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, off+d, &rdof));
469cb1e1211SMatthew G Knepley       if (ldof > 0) {
470cb1e1211SMatthew G Knepley         /* We do not own this point */
471cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
4729566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetDof(sectionAdj, goff+d, rdof));
473cb1e1211SMatthew G Knepley       } else {
474cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
475cb1e1211SMatthew G Knepley       }
476cb1e1211SMatthew G Knepley     }
477cb1e1211SMatthew G Knepley     if (found) continue;
4789566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
4799566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
4809566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
481cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
482cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
483cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof, noff;
484cb1e1211SMatthew G Knepley 
485cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
4869566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
4879566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
4889566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, padj, &noff));
489cb1e1211SMatthew G Knepley       for (d = goff; d < goff+dof-cdof; ++d) {
4909566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(sectionAdj, d, ndof-ncdof));
491cb1e1211SMatthew G Knepley       }
492cb1e1211SMatthew G Knepley     }
4939566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
49476185916SToby Isaac     if (anDof) {
49576185916SToby Isaac       for (d = goff; d < goff+dof-cdof; ++d) {
4969566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(sectionAdj, d, anDof));
49776185916SToby Isaac       }
49876185916SToby Isaac     }
499cb1e1211SMatthew G Knepley   }
5009566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(sectionAdj));
501cb1e1211SMatthew G Knepley   if (debug) {
5029566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation:\n"));
5039566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(sectionAdj, NULL));
504cb1e1211SMatthew G Knepley   }
505cb1e1211SMatthew G Knepley   /* Get adjacent indices */
5069566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(sectionAdj, &numCols));
5079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numCols, &cols));
508cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
50976185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff;
510cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
511cb1e1211SMatthew G Knepley 
5129566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
5139566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
5149566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
5159566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
516cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
517cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
518cb1e1211SMatthew G Knepley 
5199566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSectionAdj, off+d, &ldof));
5209566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, off+d, &rdof));
521cb1e1211SMatthew G Knepley       if (ldof > 0) {
522cb1e1211SMatthew G Knepley         /* We do not own this point */
523cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
524cb1e1211SMatthew G Knepley         PetscInt aoff, roff;
525cb1e1211SMatthew G Knepley 
5269566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, goff+d, &aoff));
5279566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(rootSectionAdj, off+d, &roff));
5289566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&cols[aoff], &rootAdj[roff], rdof));
529cb1e1211SMatthew G Knepley       } else {
530cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
531cb1e1211SMatthew G Knepley       }
532cb1e1211SMatthew G Knepley     }
533cb1e1211SMatthew G Knepley     if (found) continue;
5349566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
5359566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
5369566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
537cb1e1211SMatthew G Knepley     for (d = goff; d < goff+dof-cdof; ++d) {
538cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i = 0;
539cb1e1211SMatthew G Knepley 
5409566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, d, &adof));
5419566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(sectionAdj, d, &aoff));
542cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
543cb1e1211SMatthew G Knepley         const PetscInt  padj = tmpAdj[q];
544cb1e1211SMatthew G Knepley         PetscInt        ndof, ncdof, ngoff, nd;
545cb1e1211SMatthew G Knepley         const PetscInt *ncind;
546cb1e1211SMatthew G Knepley 
547cb1e1211SMatthew G Knepley         /* Adjacent points may not be in the section chart */
548cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
5499566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
5509566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
5519566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintIndices(section, padj, &ncind));
5529566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
553cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd, ++i) {
554cb1e1211SMatthew G Knepley           cols[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
555cb1e1211SMatthew G Knepley         }
556cb1e1211SMatthew G Knepley       }
55776185916SToby Isaac       for (q = 0; q < anDof; q++, i++) {
55876185916SToby Isaac         cols[aoff+i] = anchorAdj[anOff + q];
55976185916SToby Isaac       }
5602c71b3e2SJacob Faibussowitsch       PetscCheckFalse(i != adof,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of entries %D != %D for dof %D (point %D)", i, adof, d, p);
561cb1e1211SMatthew G Knepley     }
562cb1e1211SMatthew G Knepley   }
5639566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&anchorSectionAdj));
5649566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&leafSectionAdj));
5659566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&rootSectionAdj));
5669566063dSJacob Faibussowitsch   PetscCall(PetscFree(anchorAdj));
5679566063dSJacob Faibussowitsch   PetscCall(PetscFree(rootAdj));
5689566063dSJacob Faibussowitsch   PetscCall(PetscFree(tmpAdj));
569cb1e1211SMatthew G Knepley   /* Debugging */
570cb1e1211SMatthew G Knepley   if (debug) {
571cb1e1211SMatthew G Knepley     IS tmp;
5729566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Column indices\n"));
5739566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp));
5749566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
5759566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
576cb1e1211SMatthew G Knepley   }
577a9fb6443SMatthew G. Knepley 
578a9fb6443SMatthew G. Knepley   *sA     = sectionAdj;
579a9fb6443SMatthew G. Knepley   *colIdx = cols;
580a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
581a9fb6443SMatthew G. Knepley }
582a9fb6443SMatthew G. Knepley 
58325e402d2SMatthew G. Knepley 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[])
584a9fb6443SMatthew G. Knepley {
585e101f074SMatthew G. Knepley   PetscSection   section;
586a9fb6443SMatthew G. Knepley   PetscInt       rStart, rEnd, r, pStart, pEnd, p;
587a9fb6443SMatthew G. Knepley 
588a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
589a9fb6443SMatthew G. Knepley   /* This loop needs to change to a loop over points, then field dofs, which means we need to look both sections */
5909566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
5912c71b3e2SJacob Faibussowitsch   PetscCheckFalse(rStart%bs || rEnd%bs,PetscObjectComm((PetscObject) rLayout), PETSC_ERR_ARG_WRONG, "Invalid layout [%d, %d) for matrix, must be divisible by block size %d", rStart, rEnd, bs);
592a9fb6443SMatthew G. Knepley   if (f >= 0 && bs == 1) {
5939566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
5949566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
595a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
596294b7009SMatthew G. Knepley       PetscInt rS, rE;
597a9fb6443SMatthew G. Knepley 
5989566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE));
599a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
600a9fb6443SMatthew G. Knepley         PetscInt numCols, cStart, c;
601a9fb6443SMatthew G. Knepley 
6029566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
6039566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
604a9fb6443SMatthew G. Knepley         for (c = cStart; c < cStart+numCols; ++c) {
605a9fb6443SMatthew G. Knepley           if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
606a9fb6443SMatthew G. Knepley             ++dnz[r-rStart];
607a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++dnzu[r-rStart];
608a9fb6443SMatthew G. Knepley           } else {
609a9fb6443SMatthew G. Knepley             ++onz[r-rStart];
610a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++onzu[r-rStart];
611a9fb6443SMatthew G. Knepley           }
612a9fb6443SMatthew G. Knepley         }
613a9fb6443SMatthew G. Knepley       }
614a9fb6443SMatthew G. Knepley     }
615a9fb6443SMatthew G. Knepley   } else {
616cb1e1211SMatthew G Knepley     /* Only loop over blocks of rows */
617cb1e1211SMatthew G Knepley     for (r = rStart/bs; r < rEnd/bs; ++r) {
618cb1e1211SMatthew G Knepley       const PetscInt row = r*bs;
619cb1e1211SMatthew G Knepley       PetscInt       numCols, cStart, c;
620cb1e1211SMatthew G Knepley 
6219566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, row, &numCols));
6229566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(sectionAdj, row, &cStart));
623cb1e1211SMatthew G Knepley       for (c = cStart; c < cStart+numCols; ++c) {
624e7bcfa36SMatthew G. Knepley         if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
625e7bcfa36SMatthew G. Knepley           ++dnz[r-rStart/bs];
626e7bcfa36SMatthew G. Knepley           if (cols[c] >= row) ++dnzu[r-rStart/bs];
627cb1e1211SMatthew G Knepley         } else {
628e7bcfa36SMatthew G. Knepley           ++onz[r-rStart/bs];
629e7bcfa36SMatthew G. Knepley           if (cols[c] >= row) ++onzu[r-rStart/bs];
630cb1e1211SMatthew G Knepley         }
631cb1e1211SMatthew G Knepley       }
632cb1e1211SMatthew G Knepley     }
633a9fb6443SMatthew G. Knepley     for (r = 0; r < (rEnd - rStart)/bs; ++r) {
634cb1e1211SMatthew G Knepley       dnz[r]  /= bs;
635cb1e1211SMatthew G Knepley       onz[r]  /= bs;
636cb1e1211SMatthew G Knepley       dnzu[r] /= bs;
637cb1e1211SMatthew G Knepley       onzu[r] /= bs;
638cb1e1211SMatthew G Knepley     }
639cb1e1211SMatthew G Knepley   }
640a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
641a9fb6443SMatthew G. Knepley }
642a9fb6443SMatthew G. Knepley 
64325e402d2SMatthew G. Knepley static PetscErrorCode DMPlexFillMatrix_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], Mat A)
644a9fb6443SMatthew G. Knepley {
645e101f074SMatthew G. Knepley   PetscSection   section;
646a9fb6443SMatthew G. Knepley   PetscScalar   *values;
647a9fb6443SMatthew G. Knepley   PetscInt       rStart, rEnd, r, pStart, pEnd, p, len, maxRowLen = 0;
648a9fb6443SMatthew G. Knepley 
649a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
6509566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
651a9fb6443SMatthew G. Knepley   for (r = rStart; r < rEnd; ++r) {
6529566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(sectionAdj, r, &len));
653a9fb6443SMatthew G. Knepley     maxRowLen = PetscMax(maxRowLen, len);
654a9fb6443SMatthew G. Knepley   }
6559566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(maxRowLen, &values));
656a9fb6443SMatthew G. Knepley   if (f >=0 && bs == 1) {
6579566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
6589566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
659a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
660294b7009SMatthew G. Knepley       PetscInt rS, rE;
661a9fb6443SMatthew G. Knepley 
6629566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE));
663a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
6647e01ee02SMatthew G. Knepley         PetscInt numCols, cStart;
665a9fb6443SMatthew G. Knepley 
6669566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
6679566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
6689566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
669a9fb6443SMatthew G. Knepley       }
670a9fb6443SMatthew G. Knepley     }
671a9fb6443SMatthew G. Knepley   } else {
672a9fb6443SMatthew G. Knepley     for (r = rStart; r < rEnd; ++r) {
673a9fb6443SMatthew G. Knepley       PetscInt numCols, cStart;
674a9fb6443SMatthew G. Knepley 
6759566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
6769566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
6779566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
678a9fb6443SMatthew G. Knepley     }
679a9fb6443SMatthew G. Knepley   }
6809566063dSJacob Faibussowitsch   PetscCall(PetscFree(values));
681a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
682a9fb6443SMatthew G. Knepley }
683a9fb6443SMatthew G. Knepley 
68425e402d2SMatthew G. Knepley /*@C
68525e402d2SMatthew G. Knepley   DMPlexPreallocateOperator - Calculate the matrix nonzero pattern based upon the information in the DM,
68625e402d2SMatthew G. Knepley   the PetscDS it contains, and the default PetscSection.
68725e402d2SMatthew G. Knepley 
68825e402d2SMatthew G. Knepley   Collective
68925e402d2SMatthew G. Knepley 
6904165533cSJose E. Roman   Input Parameters:
69125e402d2SMatthew G. Knepley + dm   - The DMPlex
69225e402d2SMatthew G. Knepley . bs   - The matrix blocksize
69325e402d2SMatthew G. Knepley . dnz  - An array to hold the number of nonzeros in the diagonal block
69425e402d2SMatthew G. Knepley . onz  - An array to hold the number of nonzeros in the off-diagonal block
69525e402d2SMatthew G. Knepley . dnzu - An array to hold the number of nonzeros in the upper triangle of the diagonal block
69625e402d2SMatthew G. Knepley . onzu - An array to hold the number of nonzeros in the upper triangle of the off-diagonal block
697a2b725a8SWilliam Gropp - fillMatrix - If PETSC_TRUE, fill the matrix with zeros
69825e402d2SMatthew G. Knepley 
6994165533cSJose E. Roman   Output Parameter:
70025e402d2SMatthew G. Knepley . A - The preallocated matrix
70125e402d2SMatthew G. Knepley 
70225e402d2SMatthew G. Knepley   Level: advanced
70325e402d2SMatthew G. Knepley 
70425e402d2SMatthew G. Knepley .seealso: DMCreateMatrix()
70525e402d2SMatthew G. Knepley @*/
706e101f074SMatthew G. Knepley PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix)
707a9fb6443SMatthew G. Knepley {
708a9fb6443SMatthew G. Knepley   MPI_Comm       comm;
709a9fb6443SMatthew G. Knepley   PetscDS        prob;
710a9fb6443SMatthew G. Knepley   MatType        mtype;
711a9fb6443SMatthew G. Knepley   PetscSF        sf, sfDof;
712e101f074SMatthew G. Knepley   PetscSection   section;
713a9fb6443SMatthew G. Knepley   PetscInt      *remoteOffsets;
714a9fb6443SMatthew G. Knepley   PetscSection   sectionAdj[4] = {NULL, NULL, NULL, NULL};
715a9fb6443SMatthew G. Knepley   PetscInt      *cols[4]       = {NULL, NULL, NULL, NULL};
716a9fb6443SMatthew G. Knepley   PetscBool      useCone, useClosure;
7177e01ee02SMatthew G. Knepley   PetscInt       Nf, f, idx, locRows;
718a9fb6443SMatthew G. Knepley   PetscLayout    rLayout;
719a9fb6443SMatthew G. Knepley   PetscBool      isSymBlock, isSymSeqBlock, isSymMPIBlock, debug = PETSC_FALSE;
7206cded2eaSMatthew G. Knepley   PetscMPIInt    size;
721a9fb6443SMatthew G. Knepley 
722a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
723a9fb6443SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
724064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(A, MAT_CLASSID, 7);
725dadcf809SJacob Faibussowitsch   if (dnz) PetscValidIntPointer(dnz,3);
726dadcf809SJacob Faibussowitsch   if (onz) PetscValidIntPointer(onz,4);
727dadcf809SJacob Faibussowitsch   if (dnzu) PetscValidIntPointer(dnzu,5);
728dadcf809SJacob Faibussowitsch   if (onzu) PetscValidIntPointer(onzu,6);
7299566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
7309566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
7319566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
7329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL, "-dm_view_preallocation", &debug, NULL));
7339566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
7349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
7359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_Preallocate,dm,0,0,0));
736a9fb6443SMatthew G. Knepley   /* Create dof SF based on point SF */
737a9fb6443SMatthew G. Knepley   if (debug) {
738e101f074SMatthew G. Knepley     PetscSection section, sectionGlobal;
739a9fb6443SMatthew G. Knepley     PetscSF      sf;
740a9fb6443SMatthew G. Knepley 
7419566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
7429566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
7439566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
7449566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Input Section for Preallocation:\n"));
7459566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(section, NULL));
7469566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Input Global Section for Preallocation:\n"));
7479566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(sectionGlobal, NULL));
7486cded2eaSMatthew G. Knepley     if (size > 1) {
7499566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Input SF for Preallocation:\n"));
7509566063dSJacob Faibussowitsch       PetscCall(PetscSFView(sf, NULL));
751a9fb6443SMatthew G. Knepley     }
7526cded2eaSMatthew G. Knepley   }
7539566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets));
7549566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof));
7559566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
7566cded2eaSMatthew G. Knepley   if (debug && size > 1) {
7579566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Dof SF for Preallocation:\n"));
7589566063dSJacob Faibussowitsch     PetscCall(PetscSFView(sfDof, NULL));
759a9fb6443SMatthew G. Knepley   }
760a9fb6443SMatthew G. Knepley   /* Create allocation vectors from adjacency graph */
7619566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &locRows, NULL));
7629566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
7639566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
7649566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
7659566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
766a9fb6443SMatthew G. Knepley   /* There are 4 types of adjacency */
7679566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
768a9fb6443SMatthew G. Knepley   if (Nf < 1 || bs > 1) {
7699566063dSJacob Faibussowitsch     PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
770a9fb6443SMatthew G. Knepley     idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7719566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]));
7729566063dSJacob Faibussowitsch     PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
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       if (!sectionAdj[idx]) PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]));
7789566063dSJacob Faibussowitsch       PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
779a9fb6443SMatthew G. Knepley     }
780a9fb6443SMatthew G. Knepley   }
7819566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfDof));
782cb1e1211SMatthew G Knepley   /* Set matrix pattern */
7839566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu));
7849566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
78589545effSMatthew G. Knepley   /* Check for symmetric storage */
7869566063dSJacob Faibussowitsch   PetscCall(MatGetType(A, &mtype));
7879566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATSBAIJ, &isSymBlock));
7889566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock));
7899566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock));
7909566063dSJacob Faibussowitsch   if (isSymBlock || isSymSeqBlock || isSymMPIBlock) PetscCall(MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
791cb1e1211SMatthew G Knepley   /* Fill matrix with zeros */
792cb1e1211SMatthew G Knepley   if (fillMatrix) {
793a9fb6443SMatthew G. Knepley     if (Nf < 1 || bs > 1) {
7949566063dSJacob Faibussowitsch       PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
795a9fb6443SMatthew G. Knepley       idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7969566063dSJacob Faibussowitsch       PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], A));
797a9fb6443SMatthew G. Knepley     } else {
798a9fb6443SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
7999566063dSJacob Faibussowitsch         PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
800a9fb6443SMatthew G. Knepley         idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
8019566063dSJacob Faibussowitsch         PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], A));
802cb1e1211SMatthew G Knepley       }
803cb1e1211SMatthew G Knepley     }
8049566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
8059566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
806cb1e1211SMatthew G Knepley   }
8079566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
8089566063dSJacob Faibussowitsch   for (idx = 0; idx < 4; ++idx) {PetscCall(PetscSectionDestroy(&sectionAdj[idx])); PetscCall(PetscFree(cols[idx]));}
8099566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_Preallocate,dm,0,0,0));
810cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
811cb1e1211SMatthew G Knepley }
812cb1e1211SMatthew G Knepley 
813cb1e1211SMatthew G Knepley #if 0
814cb1e1211SMatthew 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)
815cb1e1211SMatthew G Knepley {
816cb1e1211SMatthew G Knepley   PetscInt       *tmpClosure,*tmpAdj,*visits;
817cb1e1211SMatthew G Knepley   PetscInt        c,cStart,cEnd,pStart,pEnd;
818cb1e1211SMatthew G Knepley 
819cb1e1211SMatthew G Knepley   PetscFunctionBegin;
8209566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
8219566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
8229566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
823cb1e1211SMatthew G Knepley 
824e7b80042SMatthew G. Knepley   maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1));
825cb1e1211SMatthew G Knepley 
8269566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
827cb1e1211SMatthew G Knepley   npoints = pEnd - pStart;
828cb1e1211SMatthew G Knepley 
8299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits));
8309566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(lvisits,pEnd-pStart));
8319566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(visits,pEnd-pStart));
8329566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
833cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
834cb1e1211SMatthew G Knepley     PetscInt *support = tmpClosure;
8359566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support));
836cb1e1211SMatthew G Knepley     for (p=0; p<supportSize; p++) lvisits[support[p]]++;
837cb1e1211SMatthew G Knepley   }
8389566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM));
8399566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd  (sf,MPIU_INT,lvisits,visits,MPI_SUM));
8409566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits,MPI_REPLACE));
8419566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd  (sf,MPIU_INT,visits,lvisits));
842cb1e1211SMatthew G Knepley 
8439566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks());
844cb1e1211SMatthew G Knepley 
8459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner));
846cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
8479566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(cellmat,maxClosureSize*maxClosureSize));
848cb1e1211SMatthew G Knepley     /*
849cb1e1211SMatthew G Knepley      Depth-first walk of transitive closure.
850cb1e1211SMatthew 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.
851cb1e1211SMatthew G Knepley      This contribution is added to dnz if owning ranks of p and q match, to onz otherwise.
852cb1e1211SMatthew G Knepley      */
853cb1e1211SMatthew G Knepley   }
854cb1e1211SMatthew G Knepley 
8559566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM));
8569566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd  (sf,MPIU_INT,lonz,onz,MPI_SUM));
857cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
858cb1e1211SMatthew G Knepley }
859cb1e1211SMatthew G Knepley #endif
860